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