xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex3.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] ="Model Equations for Advection-Diffusion\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*
5*c4762a1bSJed Brown     Page 9, Section 1.2 Model Equations for Advection-Diffusion
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown           u_t = a u_x + d u_xx
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown    The initial conditions used here different then in the book.
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown */
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown /*
14*c4762a1bSJed Brown      Helpful runtime linear solver options:
15*c4762a1bSJed Brown            -pc_type mg -da_refine 2 -snes_monitor -ksp_monitor -ts_view   (geometric multigrid with three levels)
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown */
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown /*
20*c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this file
21*c4762a1bSJed Brown    automatically includes:
22*c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h  - vectors
23*c4762a1bSJed Brown      petscmat.h  - matrices
24*c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
25*c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h   - preconditioners
26*c4762a1bSJed Brown      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
27*c4762a1bSJed Brown */
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown #include <petscts.h>
30*c4762a1bSJed Brown #include <petscdm.h>
31*c4762a1bSJed Brown #include <petscdmda.h>
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown /*
34*c4762a1bSJed Brown    User-defined application context - contains data needed by the
35*c4762a1bSJed Brown    application-provided call-back routines.
36*c4762a1bSJed Brown */
37*c4762a1bSJed Brown typedef struct {
38*c4762a1bSJed Brown   PetscScalar a,d;   /* advection and diffusion strength */
39*c4762a1bSJed Brown   PetscBool   upwind;
40*c4762a1bSJed Brown } AppCtx;
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown /*
43*c4762a1bSJed Brown    User-defined routines
44*c4762a1bSJed Brown */
45*c4762a1bSJed Brown extern PetscErrorCode InitialConditions(TS,Vec,AppCtx*);
46*c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
47*c4762a1bSJed Brown extern PetscErrorCode Solution(TS,PetscReal,Vec,AppCtx*);
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown int main(int argc,char **argv)
50*c4762a1bSJed Brown {
51*c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
52*c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
53*c4762a1bSJed Brown   Vec            U;                      /* approximate solution vector */
54*c4762a1bSJed Brown   PetscErrorCode ierr;
55*c4762a1bSJed Brown   PetscReal      dt;
56*c4762a1bSJed Brown   DM             da;
57*c4762a1bSJed Brown   PetscInt       M;
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60*c4762a1bSJed Brown      Initialize program and set problem parameters
61*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown   ierr          = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
64*c4762a1bSJed Brown   appctx.a      = 1.0;
65*c4762a1bSJed Brown   appctx.d      = 0.0;
66*c4762a1bSJed Brown   ierr          = PetscOptionsGetScalar(NULL,NULL,"-a",&appctx.a,NULL);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr          = PetscOptionsGetScalar(NULL,NULL,"-d",&appctx.d,NULL);CHKERRQ(ierr);
68*c4762a1bSJed Brown   appctx.upwind = PETSC_TRUE;
69*c4762a1bSJed Brown   ierr          = PetscOptionsGetBool(NULL,NULL,"-upwind",&appctx.upwind,NULL);CHKERRQ(ierr);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC, 60, 1, 1,NULL,&da);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
74*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75*c4762a1bSJed Brown      Create vector data structures
76*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   /*
79*c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
80*c4762a1bSJed Brown   */
81*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&U);CHKERRQ(ierr);
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84*c4762a1bSJed Brown      Create timestepping solver context
85*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   /*
91*c4762a1bSJed Brown       For linear problems with a time-dependent f(U,t) in the equation
92*c4762a1bSJed Brown      u_t = f(u,t), the user provides the discretized right-hand-side
93*c4762a1bSJed Brown       as a time-dependent matrix.
94*c4762a1bSJed Brown   */
95*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSMatrixHeat,&appctx);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = TSSetSolutionFunction(ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void*))Solution,&appctx);CHKERRQ(ierr);
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100*c4762a1bSJed Brown      Customize timestepping solver:
101*c4762a1bSJed Brown        - Set timestepping duration info
102*c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
103*c4762a1bSJed Brown      For example,
104*c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
105*c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
106*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
109*c4762a1bSJed Brown   dt   = .48/(M*M);
110*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,1000);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,100.0);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
116*c4762a1bSJed Brown 
117*c4762a1bSJed Brown   /*
118*c4762a1bSJed Brown      Evaluate initial conditions
119*c4762a1bSJed Brown   */
120*c4762a1bSJed Brown   ierr = InitialConditions(ts,U,&appctx);CHKERRQ(ierr);
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown   /*
123*c4762a1bSJed Brown      Run the timestepping solver
124*c4762a1bSJed Brown   */
125*c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
130*c4762a1bSJed Brown      are no longer needed.
131*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown   /*
138*c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
139*c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
140*c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
141*c4762a1bSJed Brown          options are chosen (e.g., -log_view).
142*c4762a1bSJed Brown   */
143*c4762a1bSJed Brown   ierr = PetscFinalize();
144*c4762a1bSJed Brown   return ierr;
145*c4762a1bSJed Brown }
146*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
147*c4762a1bSJed Brown /*
148*c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown    Input Parameter:
151*c4762a1bSJed Brown    u - uninitialized solution vector (global)
152*c4762a1bSJed Brown    appctx - user-defined application context
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown    Output Parameter:
155*c4762a1bSJed Brown    u - vector with solution at initial time (global)
156*c4762a1bSJed Brown */
157*c4762a1bSJed Brown PetscErrorCode InitialConditions(TS ts,Vec U,AppCtx *appctx)
158*c4762a1bSJed Brown {
159*c4762a1bSJed Brown   PetscScalar    *u,h;
160*c4762a1bSJed Brown   PetscErrorCode ierr;
161*c4762a1bSJed Brown   PetscInt       i,mstart,mend,xm,M;
162*c4762a1bSJed Brown   DM             da;
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
167*c4762a1bSJed Brown   h    = 1.0/M;
168*c4762a1bSJed Brown   mend = mstart + xm;
169*c4762a1bSJed Brown   /*
170*c4762a1bSJed Brown     Get a pointer to vector data.
171*c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
172*c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
173*c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
174*c4762a1bSJed Brown       the array.
175*c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
176*c4762a1bSJed Brown       C version.  See the users manual for details.
177*c4762a1bSJed Brown   */
178*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   /*
181*c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
182*c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
183*c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
184*c4762a1bSJed Brown   */
185*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) u[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown   /*
188*c4762a1bSJed Brown      Restore vector
189*c4762a1bSJed Brown   */
190*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
191*c4762a1bSJed Brown   return 0;
192*c4762a1bSJed Brown }
193*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
194*c4762a1bSJed Brown /*
195*c4762a1bSJed Brown    Solution - Computes the exact solution at a given time.
196*c4762a1bSJed Brown 
197*c4762a1bSJed Brown    Input Parameters:
198*c4762a1bSJed Brown    t - current time
199*c4762a1bSJed Brown    solution - vector in which exact solution will be computed
200*c4762a1bSJed Brown    appctx - user-defined application context
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown    Output Parameter:
203*c4762a1bSJed Brown    solution - vector with the newly computed exact solution
204*c4762a1bSJed Brown */
205*c4762a1bSJed Brown PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *appctx)
206*c4762a1bSJed Brown {
207*c4762a1bSJed Brown   PetscScalar    *u,ex1,ex2,sc1,sc2,h;
208*c4762a1bSJed Brown   PetscErrorCode ierr;
209*c4762a1bSJed Brown   PetscInt       i,mstart,mend,xm,M;
210*c4762a1bSJed Brown   DM             da;
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
213*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr);
214*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
215*c4762a1bSJed Brown   h    = 1.0/M;
216*c4762a1bSJed Brown   mend = mstart + xm;
217*c4762a1bSJed Brown   /*
218*c4762a1bSJed Brown      Get a pointer to vector data.
219*c4762a1bSJed Brown   */
220*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
221*c4762a1bSJed Brown 
222*c4762a1bSJed Brown   /*
223*c4762a1bSJed Brown      Simply write the solution directly into the array locations.
224*c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
225*c4762a1bSJed Brown   */
226*c4762a1bSJed Brown   ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*appctx->d*t);
227*c4762a1bSJed Brown   ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*appctx->d*t);
228*c4762a1bSJed Brown   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
229*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) u[i] = PetscSinScalar(sc1*(PetscReal)i + appctx->a*PETSC_PI*6.*t)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i + appctx->a*PETSC_PI*2.*t)*ex2;
230*c4762a1bSJed Brown 
231*c4762a1bSJed Brown   /*
232*c4762a1bSJed Brown      Restore vector
233*c4762a1bSJed Brown   */
234*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
235*c4762a1bSJed Brown   return 0;
236*c4762a1bSJed Brown }
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
239*c4762a1bSJed Brown /*
240*c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
241*c4762a1bSJed Brown    matrix for the heat equation.
242*c4762a1bSJed Brown 
243*c4762a1bSJed Brown    Input Parameters:
244*c4762a1bSJed Brown    ts - the TS context
245*c4762a1bSJed Brown    t - current time
246*c4762a1bSJed Brown    global_in - global input vector
247*c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown    Output Parameters:
250*c4762a1bSJed Brown    AA - Jacobian matrix
251*c4762a1bSJed Brown    BB - optionally different preconditioning matrix
252*c4762a1bSJed Brown    str - flag indicating matrix structure
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown    Notes:
255*c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
256*c4762a1bSJed Brown    in Fortran as well as in C.
257*c4762a1bSJed Brown */
258*c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec U,Mat AA,Mat BB,void *ctx)
259*c4762a1bSJed Brown {
260*c4762a1bSJed Brown   Mat            A       = AA;                /* Jacobian matrix */
261*c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
262*c4762a1bSJed Brown   PetscInt       mstart, mend;
263*c4762a1bSJed Brown   PetscErrorCode ierr;
264*c4762a1bSJed Brown   PetscInt       i,idx[3],M,xm;
265*c4762a1bSJed Brown   PetscScalar    v[3],h;
266*c4762a1bSJed Brown   DM             da;
267*c4762a1bSJed Brown 
268*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
270*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr);
271*c4762a1bSJed Brown   h    = 1.0/M;
272*c4762a1bSJed Brown   mend = mstart + xm;
273*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
274*c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
275*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
276*c4762a1bSJed Brown   /*
277*c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
278*c4762a1bSJed Brown   */
279*c4762a1bSJed Brown 
280*c4762a1bSJed Brown   /* diffusion */
281*c4762a1bSJed Brown   v[0] = appctx->d/(h*h);
282*c4762a1bSJed Brown   v[1] = -2.0*appctx->d/(h*h);
283*c4762a1bSJed Brown   v[2] = appctx->d/(h*h);
284*c4762a1bSJed Brown   if (!mstart) {
285*c4762a1bSJed Brown     idx[0] = M-1; idx[1] = 0; idx[2] = 1;
286*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&mstart,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
287*c4762a1bSJed Brown     mstart++;
288*c4762a1bSJed Brown   }
289*c4762a1bSJed Brown 
290*c4762a1bSJed Brown   if (mend == M) {
291*c4762a1bSJed Brown     mend--;
292*c4762a1bSJed Brown     idx[0] = M-2; idx[1] = M-1; idx[2] = 0;
293*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&mend,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
294*c4762a1bSJed Brown   }
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown   /*
297*c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
298*c4762a1bSJed Brown      matrix one row at a time.
299*c4762a1bSJed Brown   */
300*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
301*c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
302*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
303*c4762a1bSJed Brown   }
304*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
305*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
306*c4762a1bSJed Brown 
307*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr);
308*c4762a1bSJed Brown   mend = mstart + xm;
309*c4762a1bSJed Brown   if (!appctx->upwind) {
310*c4762a1bSJed Brown     /* advection -- centered differencing */
311*c4762a1bSJed Brown     v[0] = -.5*appctx->a/(h);
312*c4762a1bSJed Brown     v[1] = .5*appctx->a/(h);
313*c4762a1bSJed Brown     if (!mstart) {
314*c4762a1bSJed Brown       idx[0] = M-1; idx[1] = 1;
315*c4762a1bSJed Brown       ierr   = MatSetValues(A,1,&mstart,2,idx,v,ADD_VALUES);CHKERRQ(ierr);
316*c4762a1bSJed Brown       mstart++;
317*c4762a1bSJed Brown     }
318*c4762a1bSJed Brown 
319*c4762a1bSJed Brown     if (mend == M) {
320*c4762a1bSJed Brown       mend--;
321*c4762a1bSJed Brown       idx[0] = M-2; idx[1] = 0;
322*c4762a1bSJed Brown       ierr   = MatSetValues(A,1,&mend,2,idx,v,ADD_VALUES);CHKERRQ(ierr);
323*c4762a1bSJed Brown     }
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown     for (i=mstart; i<mend; i++) {
326*c4762a1bSJed Brown       idx[0] = i-1; idx[1] = i+1;
327*c4762a1bSJed Brown       ierr   = MatSetValues(A,1,&i,2,idx,v,ADD_VALUES);CHKERRQ(ierr);
328*c4762a1bSJed Brown     }
329*c4762a1bSJed Brown   } else {
330*c4762a1bSJed Brown     /* advection -- upwinding */
331*c4762a1bSJed Brown     v[0] = -appctx->a/(h);
332*c4762a1bSJed Brown     v[1] = appctx->a/(h);
333*c4762a1bSJed Brown     if (!mstart) {
334*c4762a1bSJed Brown       idx[0] = 0; idx[1] = 1;
335*c4762a1bSJed Brown       ierr   = MatSetValues(A,1,&mstart,2,idx,v,ADD_VALUES);CHKERRQ(ierr);
336*c4762a1bSJed Brown       mstart++;
337*c4762a1bSJed Brown     }
338*c4762a1bSJed Brown 
339*c4762a1bSJed Brown     if (mend == M) {
340*c4762a1bSJed Brown       mend--;
341*c4762a1bSJed Brown       idx[0] = M-1; idx[1] = 0;
342*c4762a1bSJed Brown       ierr   = MatSetValues(A,1,&mend,2,idx,v,ADD_VALUES);CHKERRQ(ierr);
343*c4762a1bSJed Brown     }
344*c4762a1bSJed Brown 
345*c4762a1bSJed Brown     for (i=mstart; i<mend; i++) {
346*c4762a1bSJed Brown       idx[0] = i; idx[1] = i+1;
347*c4762a1bSJed Brown       ierr   = MatSetValues(A,1,&i,2,idx,v,ADD_VALUES);CHKERRQ(ierr);
348*c4762a1bSJed Brown     }
349*c4762a1bSJed Brown   }
350*c4762a1bSJed Brown 
351*c4762a1bSJed Brown 
352*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
353*c4762a1bSJed Brown      Complete the matrix assembly process and set some options
354*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
355*c4762a1bSJed Brown   /*
356*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
357*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
358*c4762a1bSJed Brown      Computations can be done while messages are in transition
359*c4762a1bSJed Brown      by placing code between these two statements.
360*c4762a1bSJed Brown   */
361*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
362*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
363*c4762a1bSJed Brown 
364*c4762a1bSJed Brown   /*
365*c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
366*c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
367*c4762a1bSJed Brown   */
368*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
369*c4762a1bSJed Brown   return 0;
370*c4762a1bSJed Brown }
371*c4762a1bSJed Brown 
372*c4762a1bSJed Brown 
373*c4762a1bSJed Brown /*TEST
374*c4762a1bSJed Brown 
375*c4762a1bSJed Brown    test:
376*c4762a1bSJed Brown       args: -pc_type mg -da_refine 2  -ts_view  -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3
377*c4762a1bSJed Brown       requires: double
378*c4762a1bSJed Brown 
379*c4762a1bSJed Brown    test:
380*c4762a1bSJed Brown      suffix: 2
381*c4762a1bSJed Brown      args:  -pc_type mg -da_refine 2  -ts_view  -ts_monitor_draw_solution -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3
382*c4762a1bSJed Brown      requires: x
383*c4762a1bSJed Brown      output_file: output/ex3_1.out
384*c4762a1bSJed Brown      requires: double
385*c4762a1bSJed Brown 
386*c4762a1bSJed Brown TEST*/
387