xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex3.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Model Equations for Advection-Diffusion\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown     Page 9, Section 1.2 Model Equations for Advection-Diffusion
5c4762a1bSJed Brown 
6c4762a1bSJed Brown           u_t = a u_x + d u_xx
7c4762a1bSJed Brown 
8c4762a1bSJed Brown    The initial conditions used here different then in the book.
9c4762a1bSJed Brown 
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown      Helpful runtime linear solver options:
14c4762a1bSJed Brown            -pc_type mg -da_refine 2 -snes_monitor -ksp_monitor -ts_view   (geometric multigrid with three levels)
15c4762a1bSJed Brown 
16c4762a1bSJed Brown */
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this file
20c4762a1bSJed Brown    automatically includes:
21c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h  - vectors
22c4762a1bSJed Brown      petscmat.h  - matrices
23c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
24c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h   - preconditioners
25c4762a1bSJed Brown      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
26c4762a1bSJed Brown */
27c4762a1bSJed Brown 
28c4762a1bSJed Brown #include <petscts.h>
29c4762a1bSJed Brown #include <petscdm.h>
30c4762a1bSJed Brown #include <petscdmda.h>
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown    User-defined application context - contains data needed by the
34c4762a1bSJed Brown    application-provided call-back routines.
35c4762a1bSJed Brown */
36c4762a1bSJed Brown typedef struct {
37c4762a1bSJed Brown   PetscScalar a, d; /* advection and diffusion strength */
38c4762a1bSJed Brown   PetscBool   upwind;
39c4762a1bSJed Brown } AppCtx;
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /*
42c4762a1bSJed Brown    User-defined routines
43c4762a1bSJed Brown */
44c4762a1bSJed Brown extern PetscErrorCode InitialConditions(TS, Vec, AppCtx *);
45c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *);
46c4762a1bSJed Brown extern PetscErrorCode Solution(TS, PetscReal, Vec, AppCtx *);
47c4762a1bSJed Brown 
main(int argc,char ** argv)48d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
49d71ae5a4SJacob Faibussowitsch {
50c4762a1bSJed Brown   AppCtx    appctx; /* user-defined application context */
51c4762a1bSJed Brown   TS        ts;     /* timestepping context */
52c4762a1bSJed Brown   Vec       U;      /* approximate solution vector */
53c4762a1bSJed Brown   PetscReal dt;
54c4762a1bSJed Brown   DM        da;
55c4762a1bSJed Brown   PetscInt  M;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Initialize program and set problem parameters
59c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60c4762a1bSJed Brown 
61327415f7SBarry Smith   PetscFunctionBeginUser;
62c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
63c4762a1bSJed Brown   appctx.a = 1.0;
64c4762a1bSJed Brown   appctx.d = 0.0;
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-a", &appctx.a, NULL));
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-d", &appctx.d, NULL));
67c4762a1bSJed Brown   appctx.upwind = PETSC_TRUE;
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-upwind", &appctx.upwind, NULL));
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da));
719566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
729566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
73c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74c4762a1bSJed Brown      Create vector data structures
75c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /*
78c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
79c4762a1bSJed Brown   */
809566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &U));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown      Create timestepping solver context
84c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
879566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /*
90c4762a1bSJed Brown       For linear problems with a time-dependent f(U,t) in the equation
91dd8e379bSPierre Jolivet      u_t = f(u,t), the user provides the discretized right-hand side
92c4762a1bSJed Brown       as a time-dependent matrix.
93c4762a1bSJed Brown   */
949566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
959566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSMatrixHeat, &appctx));
969566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode (*)(TS, PetscReal, Vec, void *))Solution, &appctx));
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99c4762a1bSJed Brown      Customize timestepping solver:
100c4762a1bSJed Brown        - Set timestepping duration info
101c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
102c4762a1bSJed Brown      For example,
103c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
104c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
105c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
108c4762a1bSJed Brown   dt = .48 / (M * M);
1099566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
1109566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 1000));
1119566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 100.0));
1129566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1139566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSARKIMEX));
1149566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /*
117c4762a1bSJed Brown      Evaluate initial conditions
118c4762a1bSJed Brown   */
1199566063dSJacob Faibussowitsch   PetscCall(InitialConditions(ts, U, &appctx));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /*
122c4762a1bSJed Brown      Run the timestepping solver
123c4762a1bSJed Brown   */
1249566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
128c4762a1bSJed Brown      are no longer needed.
129c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130c4762a1bSJed Brown 
1319566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /*
136c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
137c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
138c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
139c4762a1bSJed Brown          options are chosen (e.g., -log_view).
140c4762a1bSJed Brown   */
1419566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
142b122ec5aSJacob Faibussowitsch   return 0;
143c4762a1bSJed Brown }
144c4762a1bSJed Brown /* --------------------------------------------------------------------- */
145c4762a1bSJed Brown /*
146c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    Input Parameter:
149c4762a1bSJed Brown    u - uninitialized solution vector (global)
150c4762a1bSJed Brown    appctx - user-defined application context
151c4762a1bSJed Brown 
152c4762a1bSJed Brown    Output Parameter:
153c4762a1bSJed Brown    u - vector with solution at initial time (global)
154c4762a1bSJed Brown */
InitialConditions(TS ts,Vec U,AppCtx * appctx)155d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx)
156d71ae5a4SJacob Faibussowitsch {
157c4762a1bSJed Brown   PetscScalar *u, h;
158c4762a1bSJed Brown   PetscInt     i, mstart, mend, xm, M;
159c4762a1bSJed Brown   DM           da;
160c4762a1bSJed Brown 
1613ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1629566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
1649566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
165c4762a1bSJed Brown   h    = 1.0 / M;
166c4762a1bSJed Brown   mend = mstart + xm;
167c4762a1bSJed Brown   /*
168c4762a1bSJed Brown     Get a pointer to vector data.
169c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
170c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
171c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
172c4762a1bSJed Brown       the array.
173c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
174c4762a1bSJed Brown       C version.  See the users manual for details.
175c4762a1bSJed Brown   */
1769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /*
179c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
180c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
181c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
182c4762a1bSJed Brown   */
183c4762a1bSJed Brown   for (i = mstart; i < mend; i++) u[i] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h);
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /*
186c4762a1bSJed Brown      Restore vector
187c4762a1bSJed Brown   */
1889566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
1893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
190c4762a1bSJed Brown }
191c4762a1bSJed Brown /* --------------------------------------------------------------------- */
192c4762a1bSJed Brown /*
193c4762a1bSJed Brown    Solution - Computes the exact solution at a given time.
194c4762a1bSJed Brown 
195c4762a1bSJed Brown    Input Parameters:
196c4762a1bSJed Brown    t - current time
197c4762a1bSJed Brown    solution - vector in which exact solution will be computed
198c4762a1bSJed Brown    appctx - user-defined application context
199c4762a1bSJed Brown 
200c4762a1bSJed Brown    Output Parameter:
201c4762a1bSJed Brown    solution - vector with the newly computed exact solution
202c4762a1bSJed Brown */
Solution(TS ts,PetscReal t,Vec U,AppCtx * appctx)203d71ae5a4SJacob Faibussowitsch PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx)
204d71ae5a4SJacob Faibussowitsch {
205c4762a1bSJed Brown   PetscScalar *u, ex1, ex2, sc1, sc2, h;
206c4762a1bSJed Brown   PetscInt     i, mstart, mend, xm, M;
207c4762a1bSJed Brown   DM           da;
208c4762a1bSJed Brown 
2093ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2109566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
2129566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
213c4762a1bSJed Brown   h    = 1.0 / M;
214c4762a1bSJed Brown   mend = mstart + xm;
215c4762a1bSJed Brown   /*
216c4762a1bSJed Brown      Get a pointer to vector data.
217c4762a1bSJed Brown   */
2189566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /*
221c4762a1bSJed Brown      Simply write the solution directly into the array locations.
222c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
223c4762a1bSJed Brown   */
224c4762a1bSJed Brown   ex1 = PetscExpScalar(-36. * PETSC_PI * PETSC_PI * appctx->d * t);
225c4762a1bSJed Brown   ex2 = PetscExpScalar(-4. * PETSC_PI * PETSC_PI * appctx->d * t);
2269371c9d4SSatish Balay   sc1 = PETSC_PI * 6. * h;
2279371c9d4SSatish Balay   sc2 = PETSC_PI * 2. * h;
228c4762a1bSJed 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;
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /*
231c4762a1bSJed Brown      Restore vector
232c4762a1bSJed Brown   */
2339566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
235c4762a1bSJed Brown }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown /* --------------------------------------------------------------------- */
238c4762a1bSJed Brown /*
239c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
240c4762a1bSJed Brown    matrix for the heat equation.
241c4762a1bSJed Brown 
242c4762a1bSJed Brown    Input Parameters:
243c4762a1bSJed Brown    ts - the TS context
244c4762a1bSJed Brown    t - current time
245c4762a1bSJed Brown    global_in - global input vector
246c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
247c4762a1bSJed Brown 
248c4762a1bSJed Brown    Output Parameters:
249c4762a1bSJed Brown    AA - Jacobian matrix
2507addb90fSBarry Smith    BB - optionally different matrix used to construct the preconditioner
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    Notes:
253c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
254c4762a1bSJed Brown    in Fortran as well as in C.
255c4762a1bSJed Brown */
RHSMatrixHeat(TS ts,PetscReal t,Vec U,Mat AA,Mat BB,PetscCtx ctx)256*2a8381b2SBarry Smith PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec U, Mat AA, Mat BB, PetscCtx ctx)
257d71ae5a4SJacob Faibussowitsch {
258c4762a1bSJed Brown   Mat         A      = AA;            /* Jacobian matrix */
259c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
260c4762a1bSJed Brown   PetscInt    mstart, mend;
261c4762a1bSJed Brown   PetscInt    i, idx[3], M, xm;
262c4762a1bSJed Brown   PetscScalar v[3], h;
263c4762a1bSJed Brown   DM          da;
264c4762a1bSJed Brown 
2653ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2669566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2679566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
2689566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
269c4762a1bSJed Brown   h    = 1.0 / M;
270c4762a1bSJed Brown   mend = mstart + xm;
271c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
273c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274c4762a1bSJed Brown   /*
275c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
276c4762a1bSJed Brown   */
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /* diffusion */
279c4762a1bSJed Brown   v[0] = appctx->d / (h * h);
280c4762a1bSJed Brown   v[1] = -2.0 * appctx->d / (h * h);
281c4762a1bSJed Brown   v[2] = appctx->d / (h * h);
282c4762a1bSJed Brown   if (!mstart) {
2839371c9d4SSatish Balay     idx[0] = M - 1;
2849371c9d4SSatish Balay     idx[1] = 0;
2859371c9d4SSatish Balay     idx[2] = 1;
2869566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &mstart, 3, idx, v, INSERT_VALUES));
287c4762a1bSJed Brown     mstart++;
288c4762a1bSJed Brown   }
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   if (mend == M) {
291c4762a1bSJed Brown     mend--;
2929371c9d4SSatish Balay     idx[0] = M - 2;
2939371c9d4SSatish Balay     idx[1] = M - 1;
2949371c9d4SSatish Balay     idx[2] = 0;
2959566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &mend, 3, idx, v, INSERT_VALUES));
296c4762a1bSJed Brown   }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   /*
299c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
300c4762a1bSJed Brown      matrix one row at a time.
301c4762a1bSJed Brown   */
302c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
3039371c9d4SSatish Balay     idx[0] = i - 1;
3049371c9d4SSatish Balay     idx[1] = i;
3059371c9d4SSatish Balay     idx[2] = i + 1;
3069566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
307c4762a1bSJed Brown   }
3089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FLUSH_ASSEMBLY));
3099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FLUSH_ASSEMBLY));
310c4762a1bSJed Brown 
3119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
312c4762a1bSJed Brown   mend = mstart + xm;
313c4762a1bSJed Brown   if (!appctx->upwind) {
314c4762a1bSJed Brown     /* advection -- centered differencing */
315c4762a1bSJed Brown     v[0] = -.5 * appctx->a / (h);
316c4762a1bSJed Brown     v[1] = .5 * appctx->a / (h);
317c4762a1bSJed Brown     if (!mstart) {
3189371c9d4SSatish Balay       idx[0] = M - 1;
3199371c9d4SSatish Balay       idx[1] = 1;
3209566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES));
321c4762a1bSJed Brown       mstart++;
322c4762a1bSJed Brown     }
323c4762a1bSJed Brown 
324c4762a1bSJed Brown     if (mend == M) {
325c4762a1bSJed Brown       mend--;
3269371c9d4SSatish Balay       idx[0] = M - 2;
3279371c9d4SSatish Balay       idx[1] = 0;
3289566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES));
329c4762a1bSJed Brown     }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown     for (i = mstart; i < mend; i++) {
3329371c9d4SSatish Balay       idx[0] = i - 1;
3339371c9d4SSatish Balay       idx[1] = i + 1;
3349566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES));
335c4762a1bSJed Brown     }
336c4762a1bSJed Brown   } else {
337c4762a1bSJed Brown     /* advection -- upwinding */
338c4762a1bSJed Brown     v[0] = -appctx->a / (h);
339c4762a1bSJed Brown     v[1] = appctx->a / (h);
340c4762a1bSJed Brown     if (!mstart) {
3419371c9d4SSatish Balay       idx[0] = 0;
3429371c9d4SSatish Balay       idx[1] = 1;
3439566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES));
344c4762a1bSJed Brown       mstart++;
345c4762a1bSJed Brown     }
346c4762a1bSJed Brown 
347c4762a1bSJed Brown     if (mend == M) {
348c4762a1bSJed Brown       mend--;
3499371c9d4SSatish Balay       idx[0] = M - 1;
3509371c9d4SSatish Balay       idx[1] = 0;
3519566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES));
352c4762a1bSJed Brown     }
353c4762a1bSJed Brown 
354c4762a1bSJed Brown     for (i = mstart; i < mend; i++) {
3559371c9d4SSatish Balay       idx[0] = i;
3569371c9d4SSatish Balay       idx[1] = i + 1;
3579566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES));
358c4762a1bSJed Brown     }
359c4762a1bSJed Brown   }
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
362c4762a1bSJed Brown      Complete the matrix assembly process and set some options
363c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
364c4762a1bSJed Brown   /*
365c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
366c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
367c4762a1bSJed Brown      Computations can be done while messages are in transition
368c4762a1bSJed Brown      by placing code between these two statements.
369c4762a1bSJed Brown   */
3709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   /*
374c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
375c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
376c4762a1bSJed Brown   */
3779566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
3783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
379c4762a1bSJed Brown }
380c4762a1bSJed Brown 
381c4762a1bSJed Brown /*TEST
382c4762a1bSJed Brown 
383c4762a1bSJed Brown    test:
384c4762a1bSJed Brown       args: -pc_type mg -da_refine 2 -ts_view -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3
385c4762a1bSJed Brown       requires: double
386c2eed0edSBarry Smith       filter: grep -v "total number of"
387c4762a1bSJed Brown 
388c4762a1bSJed Brown    test:
389c4762a1bSJed Brown       suffix: 2
390c4762a1bSJed 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
391c4762a1bSJed Brown       requires: x
392c4762a1bSJed Brown       output_file: output/ex3_1.out
393c4762a1bSJed Brown       requires: double
394c2eed0edSBarry Smith       filter: grep -v "total number of"
395c4762a1bSJed Brown 
396c4762a1bSJed Brown TEST*/
397