xref: /petsc/src/ts/tutorials/ex21.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Solves a time-dependent nonlinear PDE with lower and upper bounds on the interior grid points. Uses implicit\n\
3c4762a1bSJed Brown timestepping.  Runtime options include:\n\
4c4762a1bSJed Brown   -M <xg>, where <xg> = number of grid points\n\
5c4762a1bSJed Brown   -debug : Activate debugging printouts\n\
6c4762a1bSJed Brown   -nox   : Deactivate x-window graphics\n\
7c4762a1bSJed Brown   -ul   : lower bound\n\
8c4762a1bSJed Brown   -uh  : upper bound\n\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
12c4762a1bSJed Brown    Concepts: TS^Variational inequality nonlinear solver
13c4762a1bSJed Brown    Processors: n
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16c4762a1bSJed Brown /* ------------------------------------------------------------------------
17c4762a1bSJed Brown 
18c4762a1bSJed Brown    This is a variation of ex2.c to solve the PDE
19c4762a1bSJed Brown 
20c4762a1bSJed Brown                u * u_xx
21c4762a1bSJed Brown          u_t = ---------
22c4762a1bSJed Brown                2*(t+1)^2
23c4762a1bSJed Brown 
24c4762a1bSJed Brown     with box constraints on the interior grid points
25c4762a1bSJed Brown     ul <= u(t,x) <= uh with x != 0,1
26c4762a1bSJed Brown     on the domain 0 <= x <= 1, with boundary conditions
27c4762a1bSJed Brown          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
28c4762a1bSJed Brown     and initial condition
29c4762a1bSJed Brown          u(0,x) = 1 + x*x.
30c4762a1bSJed Brown 
31c4762a1bSJed Brown     The exact solution is:
32c4762a1bSJed Brown          u(t,x) = (1 + x*x) * (1 + t)
33c4762a1bSJed Brown 
34c4762a1bSJed Brown     We use by default the backward Euler method.
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   ------------------------------------------------------------------------- */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown /*
39c4762a1bSJed Brown    Include "petscts.h" to use the PETSc timestepping routines. Note that
40c4762a1bSJed Brown    this file automatically includes "petscsys.h" and other lower-level
41c4762a1bSJed Brown    PETSc include files.
42c4762a1bSJed Brown 
43c4762a1bSJed Brown    Include the "petscdmda.h" to allow us to use the distributed array data
44c4762a1bSJed Brown    structures to manage the parallel grid.
45c4762a1bSJed Brown */
46c4762a1bSJed Brown #include <petscts.h>
47c4762a1bSJed Brown #include <petscdm.h>
48c4762a1bSJed Brown #include <petscdmda.h>
49c4762a1bSJed Brown #include <petscdraw.h>
50c4762a1bSJed Brown 
51c4762a1bSJed Brown /*
52c4762a1bSJed Brown    User-defined application context - contains data needed by the
53c4762a1bSJed Brown    application-provided callback routines.
54c4762a1bSJed Brown */
55c4762a1bSJed Brown typedef struct {
56c4762a1bSJed Brown   MPI_Comm  comm;           /* communicator */
57c4762a1bSJed Brown   DM        da;             /* distributed array data structure */
58c4762a1bSJed Brown   Vec       localwork;      /* local ghosted work vector */
59c4762a1bSJed Brown   Vec       u_local;        /* local ghosted approximate solution vector */
60c4762a1bSJed Brown   Vec       solution;       /* global exact solution vector */
61c4762a1bSJed Brown   PetscInt  m;              /* total number of grid points */
62c4762a1bSJed Brown   PetscReal h;              /* mesh width: h = 1/(m-1) */
63c4762a1bSJed Brown   PetscBool debug;          /* flag (1 indicates activation of debugging printouts) */
64c4762a1bSJed Brown } AppCtx;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /*
67c4762a1bSJed Brown    User-defined routines, provided below.
68c4762a1bSJed Brown */
69c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
70c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
71c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
72c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
73c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
74c4762a1bSJed Brown extern PetscErrorCode SetBounds(Vec,Vec,PetscScalar,PetscScalar,AppCtx*);
75c4762a1bSJed Brown 
76c4762a1bSJed Brown int main(int argc,char **argv)
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
79c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
80c4762a1bSJed Brown   Mat            A;                      /* Jacobian matrix data structure */
81c4762a1bSJed Brown   Vec            u;                      /* approximate solution vector */
82c4762a1bSJed Brown   Vec            r;                      /* residual vector */
83c4762a1bSJed Brown   PetscInt       time_steps_max = 1000;  /* default max timesteps */
84c4762a1bSJed Brown   PetscErrorCode ierr;
85c4762a1bSJed Brown   PetscReal      dt;
86c4762a1bSJed Brown   PetscReal      time_total_max = 100.0; /* default max total time */
87c4762a1bSJed Brown   Vec            xl,xu; /* Lower and upper bounds on variables */
88c4762a1bSJed Brown   PetscScalar    ul=0.0,uh = 3.0;
89c4762a1bSJed Brown   PetscBool      mymonitor;
90c4762a1bSJed Brown   PetscReal      bounds[] = {1.0, 3.3};
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93c4762a1bSJed Brown      Initialize program and set problem parameters
94c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
97c4762a1bSJed Brown   ierr = PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,bounds);CHKERRQ(ierr);
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
100c4762a1bSJed Brown   appctx.m    = 60;
101c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&appctx.m,NULL);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-ul",&ul,NULL);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-uh",&uh,NULL);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);CHKERRQ(ierr);
105c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);CHKERRQ(ierr);
106c4762a1bSJed Brown   appctx.h    = 1.0/(appctx.m-1.0);
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109c4762a1bSJed Brown      Create vector data structures
110c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /*
113c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
114c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
115c4762a1bSJed Brown      total grid values spread equally among all the processors.
116c4762a1bSJed Brown   */
117c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = DMSetFromOptions(appctx.da);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = DMSetUp(appctx.da);CHKERRQ(ierr);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /*
122c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
123c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
124c4762a1bSJed Brown      have the same types.
125c4762a1bSJed Brown   */
126c4762a1bSJed Brown   ierr = DMCreateGlobalVector(appctx.da,&u);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = DMCreateLocalVector(appctx.da,&appctx.u_local);CHKERRQ(ierr);
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /*
130c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
131c4762a1bSJed Brown      create global work vector for storing exact solution.
132c4762a1bSJed Brown   */
133c4762a1bSJed Brown   ierr = VecDuplicate(appctx.u_local,&appctx.localwork);CHKERRQ(ierr);
134c4762a1bSJed Brown   ierr = VecDuplicate(u,&appctx.solution);CHKERRQ(ierr);
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* Create residual vector */
137c4762a1bSJed Brown   ierr = VecDuplicate(u,&r);CHKERRQ(ierr);
138c4762a1bSJed Brown   /* Create lower and upper bound vectors */
139c4762a1bSJed Brown   ierr = VecDuplicate(u,&xl);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = VecDuplicate(u,&xu);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = SetBounds(xl,xu,ul,uh,&appctx);CHKERRQ(ierr);
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown      Create timestepping solver context; set callback routine for
145c4762a1bSJed Brown      right-hand-side function evaluation.
146c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,r,RHSFunction,&appctx);CHKERRQ(ierr);
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown      Set optional user-defined monitoring routine
154c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   if (mymonitor) {
157c4762a1bSJed Brown     ierr = TSMonitorSet(ts,Monitor,&appctx,NULL);CHKERRQ(ierr);
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161c4762a1bSJed Brown      For nonlinear problems, the user can provide a Jacobian evaluation
162c4762a1bSJed Brown      routine (or use a finite differencing approximation).
163c4762a1bSJed Brown 
164c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine.
165c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
170c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
171c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx);CHKERRQ(ierr);
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174c4762a1bSJed Brown      Set solution vector and initial timestep
175c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   dt   = appctx.h/2.0;
178c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown      Customize timestepping solver:
182c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
183c4762a1bSJed Brown        - Set timestepping duration info
184c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
185c4762a1bSJed Brown      For example,
186c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
187c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
188c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
191c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,time_steps_max);CHKERRQ(ierr);
192c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,time_total_max);CHKERRQ(ierr);
193c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
194c4762a1bSJed Brown   /* Set lower and upper bound on the solution vector for each time step */
195c4762a1bSJed Brown   ierr = TSVISetVariableBounds(ts,xl,xu);CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199c4762a1bSJed Brown      Solve the problem
200c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /*
203c4762a1bSJed Brown      Evaluate initial conditions
204c4762a1bSJed Brown   */
205c4762a1bSJed Brown   ierr = InitialConditions(u,&appctx);CHKERRQ(ierr);
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /*
208c4762a1bSJed Brown      Run the timestepping solver
209c4762a1bSJed Brown   */
210c4762a1bSJed Brown   ierr = TSSolve(ts,u);CHKERRQ(ierr);
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
214c4762a1bSJed Brown      are no longer needed.
215c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
218c4762a1bSJed Brown   ierr = VecDestroy(&xl);CHKERRQ(ierr);
219c4762a1bSJed Brown   ierr = VecDestroy(&xu);CHKERRQ(ierr);
220c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
221c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
222c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
223c4762a1bSJed Brown   ierr = DMDestroy(&appctx.da);CHKERRQ(ierr);
224c4762a1bSJed Brown   ierr = VecDestroy(&appctx.localwork);CHKERRQ(ierr);
225c4762a1bSJed Brown   ierr = VecDestroy(&appctx.solution);CHKERRQ(ierr);
226c4762a1bSJed Brown   ierr = VecDestroy(&appctx.u_local);CHKERRQ(ierr);
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /*
229c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
230c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
231c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
232c4762a1bSJed Brown          options are chosen (e.g., -log_view).
233c4762a1bSJed Brown   */
234c4762a1bSJed Brown   ierr = PetscFinalize();
235c4762a1bSJed Brown   return ierr;
236c4762a1bSJed Brown }
237c4762a1bSJed Brown /* --------------------------------------------------------------------- */
238c4762a1bSJed Brown /*
239c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
240c4762a1bSJed Brown 
241c4762a1bSJed Brown    Input Parameters:
242c4762a1bSJed Brown    u - uninitialized solution vector (global)
243c4762a1bSJed Brown    appctx - user-defined application context
244c4762a1bSJed Brown 
245c4762a1bSJed Brown    Output Parameter:
246c4762a1bSJed Brown    u - vector with solution at initial time (global)
247c4762a1bSJed Brown */
248c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
249c4762a1bSJed Brown {
250c4762a1bSJed Brown   PetscScalar    *u_localptr,h = appctx->h,x;
251c4762a1bSJed Brown   PetscInt       i,mybase,myend;
252c4762a1bSJed Brown   PetscErrorCode ierr;
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /*
255c4762a1bSJed Brown      Determine starting point of each processor's range of
256c4762a1bSJed Brown      grid values.
257c4762a1bSJed Brown   */
258c4762a1bSJed Brown   ierr = VecGetOwnershipRange(u,&mybase,&myend);CHKERRQ(ierr);
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /*
261c4762a1bSJed Brown     Get a pointer to vector data.
262c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
263c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
264c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
265c4762a1bSJed Brown       the array.
266c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
267c4762a1bSJed Brown       C version.  See the users manual for details.
268c4762a1bSJed Brown   */
269c4762a1bSJed Brown   ierr = VecGetArray(u,&u_localptr);CHKERRQ(ierr);
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /*
272c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
273c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
274c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
275c4762a1bSJed Brown   */
276c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
277c4762a1bSJed Brown     x = h*(PetscReal)i; /* current location in global grid */
278c4762a1bSJed Brown     u_localptr[i-mybase] = 1.0 + x*x;
279c4762a1bSJed Brown   }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   /*
282c4762a1bSJed Brown      Restore vector
283c4762a1bSJed Brown   */
284c4762a1bSJed Brown   ierr = VecRestoreArray(u,&u_localptr);CHKERRQ(ierr);
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /*
287c4762a1bSJed Brown      Print debugging information if desired
288c4762a1bSJed Brown   */
289c4762a1bSJed Brown   if (appctx->debug) {
290c4762a1bSJed Brown      ierr = PetscPrintf(appctx->comm,"initial guess vector\n");CHKERRQ(ierr);
291c4762a1bSJed Brown      ierr = VecView(u,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
292c4762a1bSJed Brown   }
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   return 0;
295c4762a1bSJed Brown }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown /* --------------------------------------------------------------------- */
298c4762a1bSJed Brown /*
299c4762a1bSJed Brown   SetBounds - Sets the lower and uper bounds on the interior points
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   Input parameters:
302c4762a1bSJed Brown   xl - vector of lower bounds
303c4762a1bSJed Brown   xu - vector of upper bounds
304c4762a1bSJed Brown   ul - constant lower bound for all variables
305c4762a1bSJed Brown   uh - constant upper bound for all variables
306c4762a1bSJed Brown   appctx - Application context
307c4762a1bSJed Brown  */
308c4762a1bSJed Brown PetscErrorCode SetBounds(Vec xl, Vec xu, PetscScalar ul, PetscScalar uh,AppCtx *appctx)
309c4762a1bSJed Brown {
310c4762a1bSJed Brown   PetscErrorCode    ierr;
311c4762a1bSJed Brown   PetscScalar       *l,*u;
312c4762a1bSJed Brown   PetscMPIInt       rank,size;
313c4762a1bSJed Brown   PetscInt          localsize;
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   PetscFunctionBeginUser;
316c4762a1bSJed Brown   ierr = VecSet(xl,ul);CHKERRQ(ierr);
317c4762a1bSJed Brown   ierr = VecSet(xu,uh);CHKERRQ(ierr);
318c4762a1bSJed Brown   ierr = VecGetLocalSize(xl,&localsize);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = VecGetArray(xl,&l);CHKERRQ(ierr);
320c4762a1bSJed Brown   ierr = VecGetArray(xu,&u);CHKERRQ(ierr);
321c4762a1bSJed Brown 
322*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(appctx->comm,&rank);CHKERRMPI(ierr);
323*ffc4695bSBarry Smith   ierr = MPI_Comm_size(appctx->comm,&size);CHKERRMPI(ierr);
324c4762a1bSJed Brown   if (!rank) {
325c4762a1bSJed Brown     l[0] = -PETSC_INFINITY;
326c4762a1bSJed Brown     u[0] =  PETSC_INFINITY;
327c4762a1bSJed Brown   }
328c4762a1bSJed Brown   if (rank == size-1) {
329c4762a1bSJed Brown     l[localsize-1] = -PETSC_INFINITY;
330c4762a1bSJed Brown     u[localsize-1] = PETSC_INFINITY;
331c4762a1bSJed Brown   }
332c4762a1bSJed Brown   ierr = VecRestoreArray(xl,&l);CHKERRQ(ierr);
333c4762a1bSJed Brown   ierr = VecRestoreArray(xu,&u);CHKERRQ(ierr);
334c4762a1bSJed Brown   PetscFunctionReturn(0);
335c4762a1bSJed Brown }
336c4762a1bSJed Brown 
337c4762a1bSJed Brown /* --------------------------------------------------------------------- */
338c4762a1bSJed Brown /*
339c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
340c4762a1bSJed Brown 
341c4762a1bSJed Brown    Input Parameters:
342c4762a1bSJed Brown    t - current time
343c4762a1bSJed Brown    solution - vector in which exact solution will be computed
344c4762a1bSJed Brown    appctx - user-defined application context
345c4762a1bSJed Brown 
346c4762a1bSJed Brown    Output Parameter:
347c4762a1bSJed Brown    solution - vector with the newly computed exact solution
348c4762a1bSJed Brown */
349c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
350c4762a1bSJed Brown {
351c4762a1bSJed Brown   PetscScalar    *s_localptr,h = appctx->h,x;
352c4762a1bSJed Brown   PetscInt       i,mybase,myend;
353c4762a1bSJed Brown   PetscErrorCode ierr;
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   /*
356c4762a1bSJed Brown      Determine starting and ending points of each processor's
357c4762a1bSJed Brown      range of grid values
358c4762a1bSJed Brown   */
359c4762a1bSJed Brown   ierr = VecGetOwnershipRange(solution,&mybase,&myend);CHKERRQ(ierr);
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /*
362c4762a1bSJed Brown      Get a pointer to vector data.
363c4762a1bSJed Brown   */
364c4762a1bSJed Brown   ierr = VecGetArray(solution,&s_localptr);CHKERRQ(ierr);
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   /*
367c4762a1bSJed Brown      Simply write the solution directly into the array locations.
368c4762a1bSJed Brown      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
369c4762a1bSJed Brown   */
370c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
371c4762a1bSJed Brown     x = h*(PetscReal)i;
372c4762a1bSJed Brown     s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x);
373c4762a1bSJed Brown   }
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   /*
376c4762a1bSJed Brown      Restore vector
377c4762a1bSJed Brown   */
378c4762a1bSJed Brown   ierr = VecRestoreArray(solution,&s_localptr);CHKERRQ(ierr);
379c4762a1bSJed Brown   return 0;
380c4762a1bSJed Brown }
381c4762a1bSJed Brown /* --------------------------------------------------------------------- */
382c4762a1bSJed Brown /*
383c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
384c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
385c4762a1bSJed Brown    error in two different norms.
386c4762a1bSJed Brown 
387c4762a1bSJed Brown    Input Parameters:
388c4762a1bSJed Brown    ts     - the timestep context
389c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
390c4762a1bSJed Brown             initial condition)
391c4762a1bSJed Brown    time   - the current time
392c4762a1bSJed Brown    u      - the solution at this timestep
393c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
394c4762a1bSJed Brown             In this case we use the application context which contains
395c4762a1bSJed Brown             information about the problem size, workspace and the exact
396c4762a1bSJed Brown             solution.
397c4762a1bSJed Brown */
398c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
399c4762a1bSJed Brown {
400c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
401c4762a1bSJed Brown   PetscErrorCode ierr;
402c4762a1bSJed Brown   PetscReal      en2,en2s,enmax;
403c4762a1bSJed Brown   PetscDraw      draw;
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   /*
406c4762a1bSJed Brown      We use the default X windows viewer
407c4762a1bSJed Brown              PETSC_VIEWER_DRAW_(appctx->comm)
408c4762a1bSJed Brown      that is associated with the current communicator. This saves
409c4762a1bSJed Brown      the effort of calling PetscViewerDrawOpen() to create the window.
410c4762a1bSJed Brown      Note that if we wished to plot several items in separate windows we
411c4762a1bSJed Brown      would create each viewer with PetscViewerDrawOpen() and store them in
412c4762a1bSJed Brown      the application context, appctx.
413c4762a1bSJed Brown 
414c4762a1bSJed Brown      PetscReal buffering makes graphics look better.
415c4762a1bSJed Brown   */
416c4762a1bSJed Brown   ierr = PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm),0,&draw);CHKERRQ(ierr);
417c4762a1bSJed Brown   ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr);
418c4762a1bSJed Brown   ierr = VecView(u,PETSC_VIEWER_DRAW_(appctx->comm));CHKERRQ(ierr);
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   /*
421c4762a1bSJed Brown      Compute the exact solution at this timestep
422c4762a1bSJed Brown   */
423c4762a1bSJed Brown   ierr = ExactSolution(time,appctx->solution,appctx);CHKERRQ(ierr);
424c4762a1bSJed Brown 
425c4762a1bSJed Brown   /*
426c4762a1bSJed Brown      Print debugging information if desired
427c4762a1bSJed Brown   */
428c4762a1bSJed Brown   if (appctx->debug) {
429c4762a1bSJed Brown     ierr = PetscPrintf(appctx->comm,"Computed solution vector\n");CHKERRQ(ierr);
430c4762a1bSJed Brown     ierr = VecView(u,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
431c4762a1bSJed Brown     ierr = PetscPrintf(appctx->comm,"Exact solution vector\n");CHKERRQ(ierr);
432c4762a1bSJed Brown     ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
433c4762a1bSJed Brown   }
434c4762a1bSJed Brown 
435c4762a1bSJed Brown   /*
436c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
437c4762a1bSJed Brown   */
438c4762a1bSJed Brown   ierr = VecAXPY(appctx->solution,-1.0,u);CHKERRQ(ierr);
439c4762a1bSJed Brown   ierr = VecNorm(appctx->solution,NORM_2,&en2);CHKERRQ(ierr);
440c4762a1bSJed Brown   en2s = PetscSqrtReal(appctx->h)*en2;  /* scale the 2-norm by the grid spacing */
441c4762a1bSJed Brown   ierr = VecNorm(appctx->solution,NORM_MAX,&enmax);CHKERRQ(ierr);
442c4762a1bSJed Brown 
443c4762a1bSJed Brown   /*
444c4762a1bSJed Brown      PetscPrintf() causes only the first processor in this
445c4762a1bSJed Brown      communicator to print the timestep information.
446c4762a1bSJed Brown   */
447c4762a1bSJed Brown   ierr = PetscPrintf(appctx->comm,"Timestep %D: time = %g,2-norm error = %g, max norm error = %g\n",step,(double)time,(double)en2s,(double)enmax);CHKERRQ(ierr);
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /*
450c4762a1bSJed Brown      Print debugging information if desired
451c4762a1bSJed Brown    */
452c4762a1bSJed Brown   /*  if (appctx->debug) {
453c4762a1bSJed Brown      ierr = PetscPrintf(appctx->comm,"Error vector\n");CHKERRQ(ierr);
454c4762a1bSJed Brown      ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
455c4762a1bSJed Brown    } */
456c4762a1bSJed Brown   return 0;
457c4762a1bSJed Brown }
458c4762a1bSJed Brown /* --------------------------------------------------------------------- */
459c4762a1bSJed Brown /*
460c4762a1bSJed Brown    RHSFunction - User-provided routine that evalues the right-hand-side
461c4762a1bSJed Brown    function of the ODE.  This routine is set in the main program by
462c4762a1bSJed Brown    calling TSSetRHSFunction().  We compute:
463c4762a1bSJed Brown           global_out = F(global_in)
464c4762a1bSJed Brown 
465c4762a1bSJed Brown    Input Parameters:
466c4762a1bSJed Brown    ts         - timesteping context
467c4762a1bSJed Brown    t          - current time
468c4762a1bSJed Brown    global_in  - vector containing the current iterate
469c4762a1bSJed Brown    ctx        - (optional) user-provided context for function evaluation.
470c4762a1bSJed Brown                 In this case we use the appctx defined above.
471c4762a1bSJed Brown 
472c4762a1bSJed Brown    Output Parameter:
473c4762a1bSJed Brown    global_out - vector containing the newly evaluated function
474c4762a1bSJed Brown */
475c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx)
476c4762a1bSJed Brown {
477c4762a1bSJed Brown   AppCtx            *appctx   = (AppCtx*) ctx;     /* user-defined application context */
478c4762a1bSJed Brown   DM                da        = appctx->da;        /* distributed array */
479c4762a1bSJed Brown   Vec               local_in  = appctx->u_local;   /* local ghosted input vector */
480c4762a1bSJed Brown   Vec               localwork = appctx->localwork; /* local ghosted work vector */
481c4762a1bSJed Brown   PetscErrorCode    ierr;
482c4762a1bSJed Brown   PetscInt          i,localsize;
483c4762a1bSJed Brown   PetscMPIInt       rank,size;
484c4762a1bSJed Brown   PetscScalar       *copyptr,sc;
485c4762a1bSJed Brown   const PetscScalar *localptr;
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
488c4762a1bSJed Brown      Get ready for local function computations
489c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
490c4762a1bSJed Brown   /*
491c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
492c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
493c4762a1bSJed Brown      By placing code between these two statements, computations can be
494c4762a1bSJed Brown      done while messages are in transition.
495c4762a1bSJed Brown   */
496c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
497c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
498c4762a1bSJed Brown 
499c4762a1bSJed Brown   /*
500c4762a1bSJed Brown       Access directly the values in our local INPUT work array
501c4762a1bSJed Brown   */
502c4762a1bSJed Brown   ierr = VecGetArrayRead(local_in,&localptr);CHKERRQ(ierr);
503c4762a1bSJed Brown 
504c4762a1bSJed Brown   /*
505c4762a1bSJed Brown       Access directly the values in our local OUTPUT work array
506c4762a1bSJed Brown   */
507c4762a1bSJed Brown   ierr = VecGetArray(localwork,&copyptr);CHKERRQ(ierr);
508c4762a1bSJed Brown 
509c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
510c4762a1bSJed Brown 
511c4762a1bSJed Brown   /*
512c4762a1bSJed Brown       Evaluate our function on the nodes owned by this processor
513c4762a1bSJed Brown   */
514c4762a1bSJed Brown   ierr = VecGetLocalSize(local_in,&localsize);CHKERRQ(ierr);
515c4762a1bSJed Brown 
516c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
517c4762a1bSJed Brown      Compute entries for the locally owned part
518c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
519c4762a1bSJed Brown 
520c4762a1bSJed Brown   /*
521c4762a1bSJed Brown      Handle boundary conditions: This is done by using the boundary condition
522c4762a1bSJed Brown         u(t,boundary) = g(t,boundary)
523c4762a1bSJed Brown      for some function g. Now take the derivative with respect to t to obtain
524c4762a1bSJed Brown         u_{t}(t,boundary) = g_{t}(t,boundary)
525c4762a1bSJed Brown 
526c4762a1bSJed Brown      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
527c4762a1bSJed Brown              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
528c4762a1bSJed Brown   */
529*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(appctx->comm,&rank);CHKERRMPI(ierr);
530*ffc4695bSBarry Smith   ierr = MPI_Comm_size(appctx->comm,&size);CHKERRMPI(ierr);
531c4762a1bSJed Brown   if (!rank) copyptr[0] = 1.0;
532c4762a1bSJed Brown   if (rank == size-1) copyptr[localsize-1] = (t < .5) ? 2.0 : 0.0;
533c4762a1bSJed Brown 
534c4762a1bSJed Brown   /*
535c4762a1bSJed Brown      Handle the interior nodes where the PDE is replace by finite
536c4762a1bSJed Brown      difference operators.
537c4762a1bSJed Brown   */
538c4762a1bSJed Brown   for (i=1; i<localsize-1; i++) copyptr[i] =  localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);
539c4762a1bSJed Brown 
540c4762a1bSJed Brown   /*
541c4762a1bSJed Brown      Restore vectors
542c4762a1bSJed Brown   */
543c4762a1bSJed Brown   ierr = VecRestoreArrayRead(local_in,&localptr);CHKERRQ(ierr);
544c4762a1bSJed Brown   ierr = VecRestoreArray(localwork,&copyptr);CHKERRQ(ierr);
545c4762a1bSJed Brown 
546c4762a1bSJed Brown   /*
547c4762a1bSJed Brown      Insert values from the local OUTPUT vector into the global
548c4762a1bSJed Brown      output vector
549c4762a1bSJed Brown   */
550c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out);CHKERRQ(ierr);
551c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out);CHKERRQ(ierr);
552c4762a1bSJed Brown 
553c4762a1bSJed Brown   /* Print debugging information if desired */
554c4762a1bSJed Brown   /*  if (appctx->debug) {
555c4762a1bSJed Brown      ierr = PetscPrintf(appctx->comm,"RHS function vector\n");CHKERRQ(ierr);
556c4762a1bSJed Brown      ierr = VecView(global_out,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
557c4762a1bSJed Brown    } */
558c4762a1bSJed Brown 
559c4762a1bSJed Brown   return 0;
560c4762a1bSJed Brown }
561c4762a1bSJed Brown /* --------------------------------------------------------------------- */
562c4762a1bSJed Brown /*
563c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
564c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
565c4762a1bSJed Brown 
566c4762a1bSJed Brown    Input Parameters:
567c4762a1bSJed Brown    ts - the TS context
568c4762a1bSJed Brown    t - current time
569c4762a1bSJed Brown    global_in - global input vector
570c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
571c4762a1bSJed Brown 
572c4762a1bSJed Brown    Output Parameters:
573c4762a1bSJed Brown    AA - Jacobian matrix
574c4762a1bSJed Brown    BB - optionally different preconditioning matrix
575c4762a1bSJed Brown    str - flag indicating matrix structure
576c4762a1bSJed Brown 
577c4762a1bSJed Brown   Notes:
578c4762a1bSJed Brown   RHSJacobian computes entries for the locally owned part of the Jacobian.
579c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
580c4762a1bSJed Brown      contiguous chunks of rows across the processors.
581c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
582c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
583c4762a1bSJed Brown      appropriate processor during matrix assembly).
584c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
585c4762a1bSJed Brown      using MatSetValues().
586c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
587c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
588c4762a1bSJed Brown      in Fortran as well as in C.
589c4762a1bSJed Brown */
590c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat B,void *ctx)
591c4762a1bSJed Brown {
592c4762a1bSJed Brown   AppCtx            *appctx  = (AppCtx*)ctx;    /* user-defined application context */
593c4762a1bSJed Brown   Vec               local_in = appctx->u_local;   /* local ghosted input vector */
594c4762a1bSJed Brown   DM                da       = appctx->da;        /* distributed array */
595c4762a1bSJed Brown   PetscScalar       v[3],sc;
596c4762a1bSJed Brown   const PetscScalar *localptr;
597c4762a1bSJed Brown   PetscErrorCode    ierr;
598c4762a1bSJed Brown   PetscInt          i,mstart,mend,mstarts,mends,idx[3],is;
599c4762a1bSJed Brown 
600c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
601c4762a1bSJed Brown      Get ready for local Jacobian computations
602c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
603c4762a1bSJed Brown   /*
604c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
605c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
606c4762a1bSJed Brown      By placing code between these two statements, computations can be
607c4762a1bSJed Brown      done while messages are in transition.
608c4762a1bSJed Brown   */
609c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
610c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
611c4762a1bSJed Brown 
612c4762a1bSJed Brown   /*
613c4762a1bSJed Brown      Get pointer to vector data
614c4762a1bSJed Brown   */
615c4762a1bSJed Brown   ierr = VecGetArrayRead(local_in,&localptr);CHKERRQ(ierr);
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   /*
618c4762a1bSJed Brown      Get starting and ending locally owned rows of the matrix
619c4762a1bSJed Brown   */
620c4762a1bSJed Brown   ierr   = MatGetOwnershipRange(B,&mstarts,&mends);CHKERRQ(ierr);
621c4762a1bSJed Brown   mstart = mstarts; mend = mends;
622c4762a1bSJed Brown 
623c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
624c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
625c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
626c4762a1bSJed Brown         contiguous chunks of rows across the processors.
627c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
628c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
629c4762a1bSJed Brown         appropriate processor during matrix assembly).
630c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
631c4762a1bSJed Brown       - We can set matrix entries either using either
632c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues().
633c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
634c4762a1bSJed Brown 
635c4762a1bSJed Brown   /*
636c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
637c4762a1bSJed Brown   */
638c4762a1bSJed Brown   if (mstart == 0) {
639c4762a1bSJed Brown     v[0] = 0.0;
640c4762a1bSJed Brown     ierr = MatSetValues(B,1,&mstart,1,&mstart,v,INSERT_VALUES);CHKERRQ(ierr);
641c4762a1bSJed Brown     mstart++;
642c4762a1bSJed Brown   }
643c4762a1bSJed Brown   if (mend == appctx->m) {
644c4762a1bSJed Brown     mend--;
645c4762a1bSJed Brown     v[0] = 0.0;
646c4762a1bSJed Brown     ierr = MatSetValues(B,1,&mend,1,&mend,v,INSERT_VALUES);CHKERRQ(ierr);
647c4762a1bSJed Brown   }
648c4762a1bSJed Brown 
649c4762a1bSJed Brown   /*
650c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
651c4762a1bSJed Brown      matrix one row at a time.
652c4762a1bSJed Brown   */
653c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
654c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
655c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
656c4762a1bSJed Brown     is     = i - mstart + 1;
657c4762a1bSJed Brown     v[0]   = sc*localptr[is];
658c4762a1bSJed Brown     v[1]   = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]);
659c4762a1bSJed Brown     v[2]   = sc*localptr[is];
660c4762a1bSJed Brown     ierr   = MatSetValues(B,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
661c4762a1bSJed Brown   }
662c4762a1bSJed Brown 
663c4762a1bSJed Brown   /*
664c4762a1bSJed Brown      Restore vector
665c4762a1bSJed Brown   */
666c4762a1bSJed Brown   ierr = VecRestoreArrayRead(local_in,&localptr);CHKERRQ(ierr);
667c4762a1bSJed Brown 
668c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
669c4762a1bSJed Brown      Complete the matrix assembly process and set some options
670c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
671c4762a1bSJed Brown   /*
672c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
673c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
674c4762a1bSJed Brown      Computations can be done while messages are in transition
675c4762a1bSJed Brown      by placing code between these two statements.
676c4762a1bSJed Brown   */
677c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
678c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
679c4762a1bSJed Brown 
680c4762a1bSJed Brown   /*
681c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
682c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
683c4762a1bSJed Brown   */
684c4762a1bSJed Brown   ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
685c4762a1bSJed Brown 
686c4762a1bSJed Brown   return 0;
687c4762a1bSJed Brown }
688c4762a1bSJed Brown 
689c4762a1bSJed Brown /*TEST
690c4762a1bSJed Brown 
691c4762a1bSJed Brown     test:
692c4762a1bSJed Brown       args: -snes_type vinewtonrsls -ts_type glee -mymonitor -ts_max_steps 10 -nox
693c4762a1bSJed Brown       requires: !single
694c4762a1bSJed Brown 
695c4762a1bSJed Brown TEST*/
696