xref: /petsc/src/ts/tutorials/ex6.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves a simple time-dependent linear PDE (the heat equation).\n\
3c4762a1bSJed Brown Input parameters include:\n\
4c4762a1bSJed Brown   -m <points>, where <points> = number of grid points\n\
5c4762a1bSJed Brown   -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
6c4762a1bSJed Brown   -debug              : Activate debugging printouts\n\
7c4762a1bSJed Brown   -nox                : Deactivate x-window graphics\n\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    This program solves the one-dimensional heat equation (also called the
12c4762a1bSJed Brown    diffusion equation),
13c4762a1bSJed Brown        u_t = u_xx,
14c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
15c4762a1bSJed Brown        u(t,0) = 0, u(t,1) = 0,
16c4762a1bSJed Brown    and the initial condition
17c4762a1bSJed Brown        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
18c4762a1bSJed Brown    This is a linear, second-order, parabolic equation.
19c4762a1bSJed Brown 
20c4762a1bSJed Brown    We discretize the right-hand side using finite differences with
21c4762a1bSJed Brown    uniform grid spacing h:
22c4762a1bSJed Brown        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
23c4762a1bSJed Brown    We then demonstrate time evolution using the various TS methods by
24c4762a1bSJed Brown    running the program via
25c4762a1bSJed Brown        ex3 -ts_type <timestepping solver>
26c4762a1bSJed Brown 
27c4762a1bSJed Brown    We compare the approximate solution with the exact solution, given by
28c4762a1bSJed Brown        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
29c4762a1bSJed Brown                       3*exp(-4*pi*pi*t) * sin(2*pi*x)
30c4762a1bSJed Brown 
31c4762a1bSJed Brown    Notes:
32c4762a1bSJed Brown    This code demonstrates the TS solver interface to two variants of
33c4762a1bSJed Brown    linear problems, u_t = f(u,t), namely
34c4762a1bSJed Brown      - time-dependent f:   f(u,t) is a function of t
35c4762a1bSJed Brown      - time-independent f: f(u,t) is simply f(u)
36c4762a1bSJed Brown 
37c4762a1bSJed Brown     The parallel version of this code is ts/tutorials/ex4.c
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   ------------------------------------------------------------------------- */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /*
42c4762a1bSJed Brown    Include "ts.h" so that we can use TS solvers.  Note that this file
43c4762a1bSJed Brown    automatically includes:
44c4762a1bSJed Brown      petscsys.h  - base PETSc routines   vec.h  - vectors
45c4762a1bSJed Brown      sys.h    - system routines       mat.h  - matrices
46c4762a1bSJed Brown      is.h     - index sets            ksp.h  - Krylov subspace methods
47c4762a1bSJed Brown      viewer.h - viewers               pc.h   - preconditioners
48c4762a1bSJed Brown      snes.h - nonlinear solvers
49c4762a1bSJed Brown */
50c4762a1bSJed Brown 
51c4762a1bSJed Brown #include <petscts.h>
52c4762a1bSJed Brown #include <petscdraw.h>
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*
55c4762a1bSJed Brown    User-defined application context - contains data needed by the
56c4762a1bSJed Brown    application-provided call-back routines.
57c4762a1bSJed Brown */
58c4762a1bSJed Brown typedef struct {
59c4762a1bSJed Brown   Vec         solution;         /* global exact solution vector */
60c4762a1bSJed Brown   PetscInt    m;                /* total number of grid points */
61c4762a1bSJed Brown   PetscReal   h;                /* mesh width h = 1/(m-1) */
62c4762a1bSJed Brown   PetscBool   debug;            /* flag (1 indicates activation of debugging printouts) */
63c4762a1bSJed Brown   PetscViewer viewer1, viewer2; /* viewers for the solution and error */
64c4762a1bSJed Brown   PetscReal   norm_2, norm_max; /* error norms */
65c4762a1bSJed Brown } AppCtx;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /*
68c4762a1bSJed Brown    User-defined routines
69c4762a1bSJed Brown */
70c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
71c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *);
72c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
73c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *);
74c4762a1bSJed Brown extern PetscErrorCode MyBCRoutine(TS, PetscReal, Vec, void *);
75c4762a1bSJed Brown 
769371c9d4SSatish Balay int main(int argc, char **argv) {
77c4762a1bSJed Brown   AppCtx      appctx;                 /* user-defined application context */
78c4762a1bSJed Brown   TS          ts;                     /* timestepping context */
79c4762a1bSJed Brown   Mat         A;                      /* matrix data structure */
80c4762a1bSJed Brown   Vec         u;                      /* approximate solution vector */
81c4762a1bSJed Brown   PetscReal   time_total_max = 100.0; /* default max total time */
82c4762a1bSJed Brown   PetscInt    time_steps_max = 100;   /* default max timesteps */
83c4762a1bSJed Brown   PetscDraw   draw;                   /* drawing context */
84c4762a1bSJed Brown   PetscInt    steps, m;
85c4762a1bSJed Brown   PetscMPIInt size;
86c4762a1bSJed Brown   PetscReal   dt;
87c4762a1bSJed Brown   PetscReal   ftime;
88c4762a1bSJed Brown   PetscBool   flg;
89c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90c4762a1bSJed Brown      Initialize program and set problem parameters
91c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92c4762a1bSJed Brown 
93327415f7SBarry Smith   PetscFunctionBeginUser;
949566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
95c4762a1bSJed Brown   MPI_Comm_size(PETSC_COMM_WORLD, &size);
963c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   m = 60;
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   appctx.m        = m;
103c4762a1bSJed Brown   appctx.h        = 1.0 / (m - 1.0);
104c4762a1bSJed Brown   appctx.norm_2   = 0.0;
105c4762a1bSJed Brown   appctx.norm_max = 0.0;
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving a linear TS problem on 1 processor\n"));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown      Create vector data structures
111c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /*
114c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
115c4762a1bSJed Brown   */
1169566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &u));
1179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
121c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 380, 400, 160, &appctx.viewer1));
1249566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer1, 0, &draw));
1259566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
1269566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 0, 400, 160, &appctx.viewer2));
1279566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer2, 0, &draw));
1289566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown      Create timestepping solver context
132c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1359566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_LINEAR));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Set optional user-defined monitoring routine
139c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown 
145c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
146c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
1499566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m));
1509566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1519566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-time_dependent_rhs", &flg));
154c4762a1bSJed Brown   if (flg) {
155c4762a1bSJed Brown     /*
156c4762a1bSJed Brown        For linear problems with a time-dependent f(u,t) in the equation
157c4762a1bSJed Brown        u_t = f(u,t), the user provides the discretized right-hand-side
158c4762a1bSJed Brown        as a time-dependent matrix.
159c4762a1bSJed Brown     */
1609566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1619566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, A, A, RHSMatrixHeat, &appctx));
162c4762a1bSJed Brown   } else {
163c4762a1bSJed Brown     /*
164c4762a1bSJed Brown        For linear problems with a time-independent f(u) in the equation
165c4762a1bSJed Brown        u_t = f(u), the user provides the discretized right-hand-side
166c4762a1bSJed Brown        as a matrix only once, and then sets a null matrix evaluation
167c4762a1bSJed Brown        routine.
168c4762a1bSJed Brown     */
1699566063dSJacob Faibussowitsch     PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx));
1709566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1719566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, &appctx));
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown      Set solution vector and initial timestep
176c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   dt = appctx.h * appctx.h / 2.0;
1799566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
1809566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, u));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183c4762a1bSJed Brown      Customize timestepping solver:
184c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
185c4762a1bSJed Brown        - Set timestepping duration info
186c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
187c4762a1bSJed Brown      For example,
188c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
189c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
190c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191c4762a1bSJed Brown 
1929566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
1939566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
1949566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1959566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown      Solve the problem
199c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /*
202c4762a1bSJed Brown      Evaluate initial conditions
203c4762a1bSJed Brown   */
2049566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /*
207c4762a1bSJed Brown      Run the timestepping solver
208c4762a1bSJed Brown   */
2099566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
2109566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2119566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214c4762a1bSJed Brown      View timestepping solver info
215c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216c4762a1bSJed Brown 
2179566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "avg. error (2 norm) = %g, avg. error (max norm) = %g\n", (double)(appctx.norm_2 / steps), (double)(appctx.norm_max / steps)));
2189566063dSJacob Faibussowitsch   PetscCall(TSView(ts, PETSC_VIEWER_STDOUT_SELF));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
222c4762a1bSJed Brown      are no longer needed.
223c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224c4762a1bSJed Brown 
2259566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2289566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer1));
2299566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer2));
2309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /*
233c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
234c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
235c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
236c4762a1bSJed Brown          options are chosen (e.g., -log_view).
237c4762a1bSJed Brown   */
2389566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
239b122ec5aSJacob Faibussowitsch   return 0;
240c4762a1bSJed Brown }
241c4762a1bSJed Brown /* --------------------------------------------------------------------- */
242c4762a1bSJed Brown /*
243c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
244c4762a1bSJed Brown 
245c4762a1bSJed Brown    Input Parameter:
246c4762a1bSJed Brown    u - uninitialized solution vector (global)
247c4762a1bSJed Brown    appctx - user-defined application context
248c4762a1bSJed Brown 
249c4762a1bSJed Brown    Output Parameter:
250c4762a1bSJed Brown    u - vector with solution at initial time (global)
251c4762a1bSJed Brown */
2529371c9d4SSatish Balay PetscErrorCode InitialConditions(Vec u, AppCtx *appctx) {
253c4762a1bSJed Brown   PetscScalar *u_localptr;
254c4762a1bSJed Brown   PetscInt     i;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /*
257c4762a1bSJed Brown     Get a pointer to vector data.
258c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
259c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
260c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
261c4762a1bSJed Brown       the array.
262c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
263c4762a1bSJed Brown       C version.  See the users manual for details.
264c4762a1bSJed Brown   */
2659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &u_localptr));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /*
268c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
269c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
270c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
271c4762a1bSJed Brown   */
272c4762a1bSJed Brown   for (i = 0; i < appctx->m; i++) u_localptr[i] = PetscSinReal(PETSC_PI * i * 6. * appctx->h) + 3. * PetscSinReal(PETSC_PI * i * 2. * appctx->h);
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /*
275c4762a1bSJed Brown      Restore vector
276c4762a1bSJed Brown   */
2779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &u_localptr));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /*
280c4762a1bSJed Brown      Print debugging information if desired
281c4762a1bSJed Brown   */
2821baa6e33SBarry Smith   if (appctx->debug) PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   return 0;
285c4762a1bSJed Brown }
286c4762a1bSJed Brown /* --------------------------------------------------------------------- */
287c4762a1bSJed Brown /*
288c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
289c4762a1bSJed Brown 
290c4762a1bSJed Brown    Input Parameters:
291c4762a1bSJed Brown    t - current time
292c4762a1bSJed Brown    solution - vector in which exact solution will be computed
293c4762a1bSJed Brown    appctx - user-defined application context
294c4762a1bSJed Brown 
295c4762a1bSJed Brown    Output Parameter:
296c4762a1bSJed Brown    solution - vector with the newly computed exact solution
297c4762a1bSJed Brown */
2989371c9d4SSatish Balay PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx) {
299c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
300c4762a1bSJed Brown   PetscInt     i;
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /*
303c4762a1bSJed Brown      Get a pointer to vector data.
304c4762a1bSJed Brown   */
3059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(solution, &s_localptr));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /*
308c4762a1bSJed Brown      Simply write the solution directly into the array locations.
309c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
310c4762a1bSJed Brown   */
3119371c9d4SSatish Balay   ex1 = PetscExpReal(-36. * PETSC_PI * PETSC_PI * t);
3129371c9d4SSatish Balay   ex2 = PetscExpReal(-4. * PETSC_PI * PETSC_PI * t);
3139371c9d4SSatish Balay   sc1 = PETSC_PI * 6. * h;
3149371c9d4SSatish Balay   sc2 = PETSC_PI * 2. * h;
315c4762a1bSJed Brown   for (i = 0; i < appctx->m; i++) s_localptr[i] = PetscSinReal(PetscRealPart(sc1) * (PetscReal)i) * ex1 + 3. * PetscSinReal(PetscRealPart(sc2) * (PetscReal)i) * ex2;
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   /*
318c4762a1bSJed Brown      Restore vector
319c4762a1bSJed Brown   */
3209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(solution, &s_localptr));
321c4762a1bSJed Brown   return 0;
322c4762a1bSJed Brown }
323c4762a1bSJed Brown /* --------------------------------------------------------------------- */
324c4762a1bSJed Brown /*
325c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
326c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
327c4762a1bSJed Brown    error in two different norms.
328c4762a1bSJed Brown 
329c4762a1bSJed Brown    This example also demonstrates changing the timestep via TSSetTimeStep().
330c4762a1bSJed Brown 
331c4762a1bSJed Brown    Input Parameters:
332c4762a1bSJed Brown    ts     - the timestep context
333c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
334c4762a1bSJed Brown              initial condition)
335c4762a1bSJed Brown    crtime  - the current time
336c4762a1bSJed Brown    u      - the solution at this timestep
337c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
338c4762a1bSJed Brown             In this case we use the application context which contains
339c4762a1bSJed Brown             information about the problem size, workspace and the exact
340c4762a1bSJed Brown             solution.
341c4762a1bSJed Brown */
3429371c9d4SSatish Balay PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) {
343c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
344c4762a1bSJed Brown   PetscReal norm_2, norm_max, dt, dttol;
345c4762a1bSJed Brown   PetscBool flg;
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   /*
348c4762a1bSJed Brown      View a graph of the current iterate
349c4762a1bSJed Brown   */
3509566063dSJacob Faibussowitsch   PetscCall(VecView(u, appctx->viewer2));
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   /*
353c4762a1bSJed Brown      Compute the exact solution
354c4762a1bSJed Brown   */
3559566063dSJacob Faibussowitsch   PetscCall(ExactSolution(crtime, appctx->solution, appctx));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /*
358c4762a1bSJed Brown      Print debugging information if desired
359c4762a1bSJed Brown   */
360c4762a1bSJed Brown   if (appctx->debug) {
3619566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Computed solution vector\n"));
3629566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
3639566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Exact solution vector\n"));
3649566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
365c4762a1bSJed Brown   }
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   /*
368c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
369c4762a1bSJed Brown   */
3709566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
3719566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2));
372c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h) * norm_2;
3739566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max));
374c4762a1bSJed Brown 
3759566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
376*48a46eb9SPierre Jolivet   if (norm_2 > 1.e-2) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Timestep %" PetscInt_FMT ": step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n", step, (double)dt, (double)crtime, (double)norm_2, (double)norm_max));
377c4762a1bSJed Brown   appctx->norm_2 += norm_2;
378c4762a1bSJed Brown   appctx->norm_max += norm_max;
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   dttol = .0001;
3819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-dttol", &dttol, &flg));
382c4762a1bSJed Brown   if (dt < dttol) {
383c4762a1bSJed Brown     dt *= .999;
3849566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
385c4762a1bSJed Brown   }
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /*
388c4762a1bSJed Brown      View a graph of the error
389c4762a1bSJed Brown   */
3909566063dSJacob Faibussowitsch   PetscCall(VecView(appctx->solution, appctx->viewer1));
391c4762a1bSJed Brown 
392c4762a1bSJed Brown   /*
393c4762a1bSJed Brown      Print debugging information if desired
394c4762a1bSJed Brown   */
395c4762a1bSJed Brown   if (appctx->debug) {
3969566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error vector\n"));
3979566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
398c4762a1bSJed Brown   }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   return 0;
401c4762a1bSJed Brown }
402c4762a1bSJed Brown /* --------------------------------------------------------------------- */
403c4762a1bSJed Brown /*
404c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
405c4762a1bSJed Brown    matrix for the heat equation.
406c4762a1bSJed Brown 
407c4762a1bSJed Brown    Input Parameters:
408c4762a1bSJed Brown    ts - the TS context
409c4762a1bSJed Brown    t - current time
410c4762a1bSJed Brown    global_in - global input vector
411c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
412c4762a1bSJed Brown 
413c4762a1bSJed Brown    Output Parameters:
414c4762a1bSJed Brown    AA - Jacobian matrix
415c4762a1bSJed Brown    BB - optionally different preconditioning matrix
416c4762a1bSJed Brown    str - flag indicating matrix structure
417c4762a1bSJed Brown 
418c4762a1bSJed Brown    Notes:
419c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
420c4762a1bSJed Brown    in Fortran as well as in C.
421c4762a1bSJed Brown */
4229371c9d4SSatish Balay PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec X, Mat AA, Mat BB, void *ctx) {
423c4762a1bSJed Brown   Mat         A      = AA;            /* Jacobian matrix */
424c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
425c4762a1bSJed Brown   PetscInt    mstart = 0;
426c4762a1bSJed Brown   PetscInt    mend   = appctx->m;
427c4762a1bSJed Brown   PetscInt    i, idx[3];
428c4762a1bSJed Brown   PetscScalar v[3], stwo = -2. / (appctx->h * appctx->h), sone = -.5 * stwo;
429c4762a1bSJed Brown 
430c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
431c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
432c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
433c4762a1bSJed Brown   /*
434c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
435c4762a1bSJed Brown   */
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   mstart = 0;
438c4762a1bSJed Brown   v[0]   = 1.0;
4399566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
440c4762a1bSJed Brown   mstart++;
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   mend--;
443c4762a1bSJed Brown   v[0] = 1.0;
4449566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mend, 1, &mend, v, INSERT_VALUES));
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   /*
447c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
448c4762a1bSJed Brown      matrix one row at a time.
449c4762a1bSJed Brown   */
4509371c9d4SSatish Balay   v[0] = sone;
4519371c9d4SSatish Balay   v[1] = stwo;
4529371c9d4SSatish Balay   v[2] = sone;
453c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
4549371c9d4SSatish Balay     idx[0] = i - 1;
4559371c9d4SSatish Balay     idx[1] = i;
4569371c9d4SSatish Balay     idx[2] = i + 1;
4579566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
458c4762a1bSJed Brown   }
459c4762a1bSJed Brown 
460c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
461c4762a1bSJed Brown      Complete the matrix assembly process and set some options
462c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
463c4762a1bSJed Brown   /*
464c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
465c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
466c4762a1bSJed Brown      Computations can be done while messages are in transition
467c4762a1bSJed Brown      by placing code between these two statements.
468c4762a1bSJed Brown   */
4699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
471c4762a1bSJed Brown 
472c4762a1bSJed Brown   /*
473c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
474c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
475c4762a1bSJed Brown   */
4769566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
477c4762a1bSJed Brown 
478c4762a1bSJed Brown   return 0;
479c4762a1bSJed Brown }
480c4762a1bSJed Brown /* --------------------------------------------------------------------- */
481c4762a1bSJed Brown /*
482c4762a1bSJed Brown    Input Parameters:
483c4762a1bSJed Brown    ts - the TS context
484c4762a1bSJed Brown    t - current time
485c4762a1bSJed Brown    f - function
486c4762a1bSJed Brown    ctx - optional user-defined context, as set by TSetBCFunction()
487c4762a1bSJed Brown  */
4889371c9d4SSatish Balay PetscErrorCode MyBCRoutine(TS ts, PetscReal t, Vec f, void *ctx) {
489c4762a1bSJed Brown   AppCtx      *appctx = (AppCtx *)ctx; /* user-defined application context */
490c4762a1bSJed Brown   PetscInt     m      = appctx->m;
491c4762a1bSJed Brown   PetscScalar *fa;
492c4762a1bSJed Brown 
4939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &fa));
494c4762a1bSJed Brown   fa[0]     = 0.0;
495c4762a1bSJed Brown   fa[m - 1] = 1.0;
4969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &fa));
4979566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "t=%g\n", (double)t));
498c4762a1bSJed Brown 
499c4762a1bSJed Brown   return 0;
500c4762a1bSJed Brown }
501c4762a1bSJed Brown 
502c4762a1bSJed Brown /*TEST
503c4762a1bSJed Brown 
504c4762a1bSJed Brown     test:
505c4762a1bSJed Brown       args: -nox -ts_max_steps 4
506c4762a1bSJed Brown 
507c4762a1bSJed Brown TEST*/
508