xref: /petsc/src/ts/tutorials/ex3.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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   -use_ifunc          : Use IFunction/IJacobian interface\n\
7c4762a1bSJed Brown   -debug              : Activate debugging printouts\n\
8c4762a1bSJed Brown   -nox                : Deactivate x-window graphics\n\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /* ------------------------------------------------------------------------
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    This program solves the one-dimensional heat equation (also called the
13c4762a1bSJed Brown    diffusion equation),
14c4762a1bSJed Brown        u_t = u_xx,
15c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
16c4762a1bSJed Brown        u(t,0) = 0, u(t,1) = 0,
17c4762a1bSJed Brown    and the initial condition
18c4762a1bSJed Brown        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
19c4762a1bSJed Brown    This is a linear, second-order, parabolic equation.
20c4762a1bSJed Brown 
21c4762a1bSJed Brown    We discretize the right-hand side using finite differences with
22c4762a1bSJed Brown    uniform grid spacing h:
23c4762a1bSJed Brown        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
24c4762a1bSJed Brown    We then demonstrate time evolution using the various TS methods by
25c4762a1bSJed Brown    running the program via
26c4762a1bSJed Brown        ex3 -ts_type <timestepping solver>
27c4762a1bSJed Brown 
28c4762a1bSJed Brown    We compare the approximate solution with the exact solution, given by
29c4762a1bSJed Brown        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
30c4762a1bSJed Brown                       3*exp(-4*pi*pi*t) * sin(2*pi*x)
31c4762a1bSJed Brown 
32c4762a1bSJed Brown    Notes:
33c4762a1bSJed Brown    This code demonstrates the TS solver interface to two variants of
34c4762a1bSJed Brown    linear problems, u_t = f(u,t), namely
35c4762a1bSJed Brown      - time-dependent f:   f(u,t) is a function of t
36c4762a1bSJed Brown      - time-independent f: f(u,t) is simply f(u)
37c4762a1bSJed Brown 
38c4762a1bSJed Brown     The parallel version of this code is ts/tutorials/ex4.c
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   ------------------------------------------------------------------------- */
41c4762a1bSJed Brown 
42c4762a1bSJed Brown /*
43c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this file
44c4762a1bSJed Brown    automatically includes:
45c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h  - vectors
46c4762a1bSJed Brown      petscmat.h  - matrices
47c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
48c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h   - preconditioners
49c4762a1bSJed Brown      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
50c4762a1bSJed Brown */
51c4762a1bSJed Brown 
52c4762a1bSJed Brown #include <petscts.h>
53c4762a1bSJed Brown #include <petscdraw.h>
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*
56c4762a1bSJed Brown    User-defined application context - contains data needed by the
57c4762a1bSJed Brown    application-provided call-back routines.
58c4762a1bSJed Brown */
59c4762a1bSJed Brown typedef struct {
60c4762a1bSJed Brown   Vec         solution;         /* global exact solution vector */
61c4762a1bSJed Brown   PetscInt    m;                /* total number of grid points */
62c4762a1bSJed Brown   PetscReal   h;                /* mesh width h = 1/(m-1) */
63c4762a1bSJed Brown   PetscBool   debug;            /* flag (1 indicates activation of debugging printouts) */
64c4762a1bSJed Brown   PetscViewer viewer1, viewer2; /* viewers for the solution and error */
65c4762a1bSJed Brown   PetscReal   norm_2, norm_max; /* error norms */
66c4762a1bSJed Brown   Mat         A;                /* RHS mat, used with IFunction interface */
67c4762a1bSJed Brown   PetscReal   oshift;           /* old shift applied, prevent to recompute the IJacobian */
68c4762a1bSJed Brown } AppCtx;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown /*
71c4762a1bSJed Brown    User-defined routines
72c4762a1bSJed Brown */
73c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
74c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *);
75c4762a1bSJed Brown extern PetscErrorCode IFunctionHeat(TS, PetscReal, Vec, Vec, Vec, void *);
76c4762a1bSJed Brown extern PetscErrorCode IJacobianHeat(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
77c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
78c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *);
79c4762a1bSJed Brown 
80d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
81d71ae5a4SJacob Faibussowitsch {
82c4762a1bSJed Brown   AppCtx      appctx;                 /* user-defined application context */
83c4762a1bSJed Brown   TS          ts;                     /* timestepping context */
84c4762a1bSJed Brown   Mat         A;                      /* matrix data structure */
85c4762a1bSJed Brown   Vec         u;                      /* approximate solution vector */
86c4762a1bSJed Brown   PetscReal   time_total_max = 100.0; /* default max total time */
87c4762a1bSJed Brown   PetscInt    time_steps_max = 100;   /* default max timesteps */
88c4762a1bSJed Brown   PetscDraw   draw;                   /* drawing context */
89c4762a1bSJed Brown   PetscInt    steps, m;
90c4762a1bSJed Brown   PetscMPIInt size;
91c4762a1bSJed Brown   PetscReal   dt;
92c4762a1bSJed Brown   PetscBool   flg, flg_string;
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95c4762a1bSJed Brown      Initialize program and set problem parameters
96c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97c4762a1bSJed Brown 
98327415f7SBarry Smith   PetscFunctionBeginUser;
999566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1013c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   m = 60;
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
106c4762a1bSJed Brown   flg_string = PETSC_FALSE;
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_string_viewer", &flg_string, NULL));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   appctx.m        = m;
110c4762a1bSJed Brown   appctx.h        = 1.0 / (m - 1.0);
111c4762a1bSJed Brown   appctx.norm_2   = 0.0;
112c4762a1bSJed Brown   appctx.norm_max = 0.0;
113c4762a1bSJed Brown 
1149566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving a linear TS problem on 1 processor\n"));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Create vector data structures
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /*
121c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
122c4762a1bSJed Brown   */
1239566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &u));
1249566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
128c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129c4762a1bSJed Brown 
1309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 380, 400, 160, &appctx.viewer1));
1319566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer1, 0, &draw));
1329566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
1339566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 0, 400, 160, &appctx.viewer2));
1349566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer2, 0, &draw));
1359566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Create timestepping solver context
139c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1429566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_LINEAR));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Set optional user-defined monitoring routine
146c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147c4762a1bSJed Brown 
14848a46eb9SPierre Jolivet   if (!flg_string) PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151c4762a1bSJed Brown 
152c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
153c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
1569566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m));
1579566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1589566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   flg = PETSC_FALSE;
1619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_ifunc", &flg, NULL));
162c4762a1bSJed Brown   if (!flg) {
163c4762a1bSJed Brown     appctx.A = NULL;
1649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-time_dependent_rhs", &flg, NULL));
165c4762a1bSJed Brown     if (flg) {
166c4762a1bSJed Brown       /*
167c4762a1bSJed Brown          For linear problems with a time-dependent f(u,t) in the equation
168c4762a1bSJed Brown          u_t = f(u,t), the user provides the discretized right-hand-side
169c4762a1bSJed Brown          as a time-dependent matrix.
170c4762a1bSJed Brown       */
1719566063dSJacob Faibussowitsch       PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1729566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(ts, A, A, RHSMatrixHeat, &appctx));
173c4762a1bSJed Brown     } else {
174c4762a1bSJed Brown       /*
175c4762a1bSJed Brown          For linear problems with a time-independent f(u) in the equation
176c4762a1bSJed Brown          u_t = f(u), the user provides the discretized right-hand-side
177c4762a1bSJed Brown          as a matrix only once, and then sets the special Jacobian evaluation
178c4762a1bSJed Brown          routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian.
179c4762a1bSJed Brown       */
1809566063dSJacob Faibussowitsch       PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx));
1819566063dSJacob Faibussowitsch       PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1829566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, &appctx));
183c4762a1bSJed Brown     }
184c4762a1bSJed Brown   } else {
185c4762a1bSJed Brown     Mat J;
186c4762a1bSJed Brown 
1879566063dSJacob Faibussowitsch     PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx));
1889566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &J));
1899566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunctionHeat, &appctx));
1909566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, J, J, IJacobianHeat, &appctx));
1919566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&J));
192c4762a1bSJed Brown 
1939566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
194c4762a1bSJed Brown     appctx.A      = A;
195c4762a1bSJed Brown     appctx.oshift = PETSC_MIN_REAL;
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown      Set solution vector and initial timestep
199c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   dt = appctx.h * appctx.h / 2.0;
2029566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205c4762a1bSJed Brown      Customize timestepping solver:
206c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
207c4762a1bSJed Brown        - Set timestepping duration info
208c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
209c4762a1bSJed Brown      For example,
210c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
211c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
212c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213c4762a1bSJed Brown 
2149566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
2159566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
2169566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2179566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220c4762a1bSJed Brown      Solve the problem
221c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /*
224c4762a1bSJed Brown      Evaluate initial conditions
225c4762a1bSJed Brown   */
2269566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /*
229c4762a1bSJed Brown      Run the timestepping solver
230c4762a1bSJed Brown   */
2319566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
2329566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235c4762a1bSJed Brown      View timestepping solver info
236c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
237c4762a1bSJed Brown 
2389566063dSJacob 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)));
239c4762a1bSJed Brown   if (!flg_string) {
2409566063dSJacob Faibussowitsch     PetscCall(TSView(ts, PETSC_VIEWER_STDOUT_SELF));
241c4762a1bSJed Brown   } else {
242c4762a1bSJed Brown     PetscViewer stringviewer;
243c4762a1bSJed Brown     char        string[512];
244c4762a1bSJed Brown     const char *outstring;
245c4762a1bSJed Brown 
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringOpen(PETSC_COMM_WORLD, string, sizeof(string), &stringviewer));
2479566063dSJacob Faibussowitsch     PetscCall(TSView(ts, stringviewer));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringGetStringRead(stringviewer, &outstring, NULL));
2493c633725SBarry Smith     PetscCheck((char *)outstring == (char *)string, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "String returned from viewer does not equal original string");
2509566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output from string viewer:%s\n", outstring));
2519566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&stringviewer));
252c4762a1bSJed Brown   }
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
255c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
256c4762a1bSJed Brown      are no longer needed.
257c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2629566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer1));
2639566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer2));
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
2659566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.A));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /*
268c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
269c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
270c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
271c4762a1bSJed Brown          options are chosen (e.g., -log_view).
272c4762a1bSJed Brown   */
2739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
274b122ec5aSJacob Faibussowitsch   return 0;
275c4762a1bSJed Brown }
276c4762a1bSJed Brown /* --------------------------------------------------------------------- */
277c4762a1bSJed Brown /*
278c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
279c4762a1bSJed Brown 
280c4762a1bSJed Brown    Input Parameter:
281c4762a1bSJed Brown    u - uninitialized solution vector (global)
282c4762a1bSJed Brown    appctx - user-defined application context
283c4762a1bSJed Brown 
284c4762a1bSJed Brown    Output Parameter:
285c4762a1bSJed Brown    u - vector with solution at initial time (global)
286c4762a1bSJed Brown */
287d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
288d71ae5a4SJacob Faibussowitsch {
289c4762a1bSJed Brown   PetscScalar *u_localptr, h = appctx->h;
290c4762a1bSJed Brown   PetscInt     i;
291c4762a1bSJed Brown 
292*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
293c4762a1bSJed Brown   /*
294c4762a1bSJed Brown     Get a pointer to vector data.
295c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
296c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
297c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
298c4762a1bSJed Brown       the array.
299c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
300c4762a1bSJed Brown       C version.  See the users manual for details.
301c4762a1bSJed Brown   */
3029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(u, &u_localptr));
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /*
305c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
306c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
307c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
308c4762a1bSJed Brown   */
309c4762a1bSJed Brown   for (i = 0; i < appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h);
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   /*
312c4762a1bSJed Brown      Restore vector
313c4762a1bSJed Brown   */
3149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(u, &u_localptr));
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   /*
317c4762a1bSJed Brown      Print debugging information if desired
318c4762a1bSJed Brown   */
319c4762a1bSJed Brown   if (appctx->debug) {
3209566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess vector\n"));
3219566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
322c4762a1bSJed Brown   }
323c4762a1bSJed Brown 
324*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
325c4762a1bSJed Brown }
326c4762a1bSJed Brown /* --------------------------------------------------------------------- */
327c4762a1bSJed Brown /*
328c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
329c4762a1bSJed Brown 
330c4762a1bSJed Brown    Input Parameters:
331c4762a1bSJed Brown    t - current time
332c4762a1bSJed Brown    solution - vector in which exact solution will be computed
333c4762a1bSJed Brown    appctx - user-defined application context
334c4762a1bSJed Brown 
335c4762a1bSJed Brown    Output Parameter:
336c4762a1bSJed Brown    solution - vector with the newly computed exact solution
337c4762a1bSJed Brown */
338d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx)
339d71ae5a4SJacob Faibussowitsch {
340c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2, tc = t;
341c4762a1bSJed Brown   PetscInt     i;
342c4762a1bSJed Brown 
343*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
344c4762a1bSJed Brown   /*
345c4762a1bSJed Brown      Get a pointer to vector data.
346c4762a1bSJed Brown   */
3479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(solution, &s_localptr));
348c4762a1bSJed Brown 
349c4762a1bSJed Brown   /*
350c4762a1bSJed Brown      Simply write the solution directly into the array locations.
351c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
352c4762a1bSJed Brown   */
353c4762a1bSJed Brown   ex1 = PetscExpScalar(-36. * PETSC_PI * PETSC_PI * tc);
354c4762a1bSJed Brown   ex2 = PetscExpScalar(-4. * PETSC_PI * PETSC_PI * tc);
3559371c9d4SSatish Balay   sc1 = PETSC_PI * 6. * h;
3569371c9d4SSatish Balay   sc2 = PETSC_PI * 2. * h;
357c4762a1bSJed Brown   for (i = 0; i < appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1 * (PetscReal)i) * ex1 + 3. * PetscSinScalar(sc2 * (PetscReal)i) * ex2;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   /*
360c4762a1bSJed Brown      Restore vector
361c4762a1bSJed Brown   */
3629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(solution, &s_localptr));
363*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
364c4762a1bSJed Brown }
365c4762a1bSJed Brown /* --------------------------------------------------------------------- */
366c4762a1bSJed Brown /*
367c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
368c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
369c4762a1bSJed Brown    error in two different norms.
370c4762a1bSJed Brown 
371c4762a1bSJed Brown    This example also demonstrates changing the timestep via TSSetTimeStep().
372c4762a1bSJed Brown 
373c4762a1bSJed Brown    Input Parameters:
374c4762a1bSJed Brown    ts     - the timestep context
375c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
376c4762a1bSJed Brown              initial condition)
377c4762a1bSJed Brown    time   - the current time
378c4762a1bSJed Brown    u      - the solution at this timestep
379c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
380c4762a1bSJed Brown             In this case we use the application context which contains
381c4762a1bSJed Brown             information about the problem size, workspace and the exact
382c4762a1bSJed Brown             solution.
383c4762a1bSJed Brown */
384d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx)
385d71ae5a4SJacob Faibussowitsch {
386c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
387c4762a1bSJed Brown   PetscReal norm_2, norm_max, dt, dttol;
388c4762a1bSJed Brown 
389*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
390c4762a1bSJed Brown   /*
391c4762a1bSJed Brown      View a graph of the current iterate
392c4762a1bSJed Brown   */
3939566063dSJacob Faibussowitsch   PetscCall(VecView(u, appctx->viewer2));
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   /*
396c4762a1bSJed Brown      Compute the exact solution
397c4762a1bSJed Brown   */
3989566063dSJacob Faibussowitsch   PetscCall(ExactSolution(time, appctx->solution, appctx));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   /*
401c4762a1bSJed Brown      Print debugging information if desired
402c4762a1bSJed Brown   */
403c4762a1bSJed Brown   if (appctx->debug) {
4049566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Computed solution vector\n"));
4059566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
4069566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Exact solution vector\n"));
4079566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
408c4762a1bSJed Brown   }
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   /*
411c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
412c4762a1bSJed Brown   */
4139566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
4149566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2));
415c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h) * norm_2;
4169566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max));
417c4762a1bSJed Brown 
4189566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
41963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep %3" PetscInt_FMT ": step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n", step, (double)dt, (double)time, (double)norm_2, (double)norm_max));
420c4762a1bSJed Brown 
421c4762a1bSJed Brown   appctx->norm_2 += norm_2;
422c4762a1bSJed Brown   appctx->norm_max += norm_max;
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   dttol = .0001;
4259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-dttol", &dttol, NULL));
426c4762a1bSJed Brown   if (dt < dttol) {
427c4762a1bSJed Brown     dt *= .999;
4289566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
429c4762a1bSJed Brown   }
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   /*
432c4762a1bSJed Brown      View a graph of the error
433c4762a1bSJed Brown   */
4349566063dSJacob Faibussowitsch   PetscCall(VecView(appctx->solution, appctx->viewer1));
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /*
437c4762a1bSJed Brown      Print debugging information if desired
438c4762a1bSJed Brown   */
439c4762a1bSJed Brown   if (appctx->debug) {
4409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error vector\n"));
4419566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
442c4762a1bSJed Brown   }
443c4762a1bSJed Brown 
444*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
445c4762a1bSJed Brown }
446c4762a1bSJed Brown /* --------------------------------------------------------------------- */
447c4762a1bSJed Brown /*
448c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
449c4762a1bSJed Brown    matrix for the heat equation.
450c4762a1bSJed Brown 
451c4762a1bSJed Brown    Input Parameters:
452c4762a1bSJed Brown    ts - the TS context
453c4762a1bSJed Brown    t - current time
454c4762a1bSJed Brown    global_in - global input vector
455c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
456c4762a1bSJed Brown 
457c4762a1bSJed Brown    Output Parameters:
458c4762a1bSJed Brown    AA - Jacobian matrix
459c4762a1bSJed Brown    BB - optionally different preconditioning matrix
460c4762a1bSJed Brown    str - flag indicating matrix structure
461c4762a1bSJed Brown 
462c4762a1bSJed Brown    Notes:
463c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
464c4762a1bSJed Brown    in Fortran as well as in C.
465c4762a1bSJed Brown */
466d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec X, Mat AA, Mat BB, void *ctx)
467d71ae5a4SJacob Faibussowitsch {
468c4762a1bSJed Brown   Mat         A      = AA;            /* Jacobian matrix */
469c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
470c4762a1bSJed Brown   PetscInt    mstart = 0;
471c4762a1bSJed Brown   PetscInt    mend   = appctx->m;
472c4762a1bSJed Brown   PetscInt    i, idx[3];
473c4762a1bSJed Brown   PetscScalar v[3], stwo = -2. / (appctx->h * appctx->h), sone = -.5 * stwo;
474c4762a1bSJed Brown 
475*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
476c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
477c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
478c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
479c4762a1bSJed Brown   /*
480c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
481c4762a1bSJed Brown   */
482c4762a1bSJed Brown 
483c4762a1bSJed Brown   mstart = 0;
484c4762a1bSJed Brown   v[0]   = 1.0;
4859566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
486c4762a1bSJed Brown   mstart++;
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   mend--;
489c4762a1bSJed Brown   v[0] = 1.0;
4909566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mend, 1, &mend, v, INSERT_VALUES));
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   /*
493c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
494c4762a1bSJed Brown      matrix one row at a time.
495c4762a1bSJed Brown   */
4969371c9d4SSatish Balay   v[0] = sone;
4979371c9d4SSatish Balay   v[1] = stwo;
4989371c9d4SSatish Balay   v[2] = sone;
499c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
5009371c9d4SSatish Balay     idx[0] = i - 1;
5019371c9d4SSatish Balay     idx[1] = i;
5029371c9d4SSatish Balay     idx[2] = i + 1;
5039566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
504c4762a1bSJed Brown   }
505c4762a1bSJed Brown 
506c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
507c4762a1bSJed Brown      Complete the matrix assembly process and set some options
508c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
509c4762a1bSJed Brown   /*
510c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
511c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
512c4762a1bSJed Brown      Computations can be done while messages are in transition
513c4762a1bSJed Brown      by placing code between these two statements.
514c4762a1bSJed Brown   */
5159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
517c4762a1bSJed Brown 
518c4762a1bSJed Brown   /*
519c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
520c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
521c4762a1bSJed Brown   */
5229566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
523c4762a1bSJed Brown 
524*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
525c4762a1bSJed Brown }
526c4762a1bSJed Brown 
527d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionHeat(TS ts, PetscReal t, Vec X, Vec Xdot, Vec r, void *ctx)
528d71ae5a4SJacob Faibussowitsch {
529c4762a1bSJed Brown   AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
530c4762a1bSJed Brown 
531*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
5329566063dSJacob Faibussowitsch   PetscCall(MatMult(appctx->A, X, r));
5339566063dSJacob Faibussowitsch   PetscCall(VecAYPX(r, -1.0, Xdot));
534*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535c4762a1bSJed Brown }
536c4762a1bSJed Brown 
537d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobianHeat(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal s, Mat A, Mat B, void *ctx)
538d71ae5a4SJacob Faibussowitsch {
539c4762a1bSJed Brown   AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
540c4762a1bSJed Brown 
541*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
542*3ba16761SJacob Faibussowitsch   if (appctx->oshift == s) PetscFunctionReturn(PETSC_SUCCESS);
5439566063dSJacob Faibussowitsch   PetscCall(MatCopy(appctx->A, A, SAME_NONZERO_PATTERN));
5449566063dSJacob Faibussowitsch   PetscCall(MatScale(A, -1));
5459566063dSJacob Faibussowitsch   PetscCall(MatShift(A, s));
5469566063dSJacob Faibussowitsch   PetscCall(MatCopy(A, B, SAME_NONZERO_PATTERN));
547c4762a1bSJed Brown   appctx->oshift = s;
548*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
549c4762a1bSJed Brown }
550c4762a1bSJed Brown 
551c4762a1bSJed Brown /*TEST
552c4762a1bSJed Brown 
553c4762a1bSJed Brown     test:
554c4762a1bSJed Brown       args: -nox -ts_type ssp -ts_dt 0.0005
555c4762a1bSJed Brown 
556c4762a1bSJed Brown     test:
557c4762a1bSJed Brown       suffix: 2
558c4762a1bSJed Brown       args: -nox -ts_type ssp -ts_dt 0.0005 -time_dependent_rhs 1
559c4762a1bSJed Brown 
560c4762a1bSJed Brown     test:
561c4762a1bSJed Brown       suffix: 3
562c4762a1bSJed Brown       args:  -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason
563c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
564c4762a1bSJed Brown       requires: !single
565c4762a1bSJed Brown 
566c4762a1bSJed Brown     test:
567c4762a1bSJed Brown       suffix: 4
568c4762a1bSJed Brown       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason
569c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
570c4762a1bSJed Brown 
571c4762a1bSJed Brown     test:
572c4762a1bSJed Brown       suffix: 5
573c4762a1bSJed Brown       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs
574c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
575c4762a1bSJed Brown 
576c4762a1bSJed Brown     test:
577c4762a1bSJed Brown       requires: !single
578c4762a1bSJed Brown       suffix: pod_guess
579c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason
580c4762a1bSJed Brown 
581c4762a1bSJed Brown     test:
582c4762a1bSJed Brown       requires: !single
583c4762a1bSJed Brown       suffix: pod_guess_Ainner
584c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -ksp_guess_pod_Ainner -pc_type none -ksp_converged_reason
585c4762a1bSJed Brown 
586c4762a1bSJed Brown     test:
587c4762a1bSJed Brown       requires: !single
588c4762a1bSJed Brown       suffix: fischer_guess
589c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason
590c4762a1bSJed Brown 
591c4762a1bSJed Brown     test:
592c4762a1bSJed Brown       requires: !single
593c4762a1bSJed Brown       suffix: fischer_guess_2
594c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 2,10 -pc_type none -ksp_converged_reason
595c4762a1bSJed Brown 
596c4762a1bSJed Brown     test:
597c4762a1bSJed Brown       requires: !single
598cbb17d71SDavid Wells       suffix: fischer_guess_3
599cbb17d71SDavid Wells       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 3,10 -pc_type none -ksp_converged_reason
600cbb17d71SDavid Wells 
601cbb17d71SDavid Wells     test:
602cbb17d71SDavid Wells       requires: !single
603c4762a1bSJed Brown       suffix: stringview
604c4762a1bSJed Brown       args: -nox -ts_type rosw -test_string_viewer
605c4762a1bSJed Brown 
606c4762a1bSJed Brown     test:
607c4762a1bSJed Brown       requires: !single
608c4762a1bSJed Brown       suffix: stringview_euler
609c4762a1bSJed Brown       args: -nox -ts_type euler -test_string_viewer
610c4762a1bSJed Brown 
611c4762a1bSJed Brown TEST*/
612