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