xref: /petsc/src/ts/tutorials/ex3.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Solves a simple time-dependent linear PDE (the heat equation).\n\
2c4762a1bSJed Brown Input parameters include:\n\
3c4762a1bSJed Brown   -m <points>, where <points> = number of grid points\n\
4c4762a1bSJed Brown   -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
5c4762a1bSJed Brown   -use_ifunc          : Use IFunction/IJacobian interface\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 "petscts.h" so that we can use TS solvers.  Note that this file
43c4762a1bSJed Brown    automatically includes:
44c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h  - vectors
45c4762a1bSJed Brown      petscmat.h  - matrices
46c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
47c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h   - preconditioners
48c4762a1bSJed Brown      petscksp.h   - linear solvers        petscsnes.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   Mat         A;                /* RHS mat, used with IFunction interface */
66c4762a1bSJed Brown   PetscReal   oshift;           /* old shift applied, prevent to recompute the IJacobian */
67c4762a1bSJed Brown } AppCtx;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown /*
70c4762a1bSJed Brown    User-defined routines
71c4762a1bSJed Brown */
72c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
73c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *);
74c4762a1bSJed Brown extern PetscErrorCode IFunctionHeat(TS, PetscReal, Vec, Vec, Vec, void *);
75c4762a1bSJed Brown extern PetscErrorCode IJacobianHeat(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
76c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
77c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *);
78c4762a1bSJed Brown 
main(int argc,char ** argv)79d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
80d71ae5a4SJacob Faibussowitsch {
81c4762a1bSJed Brown   AppCtx      appctx;                 /* user-defined application context */
82c4762a1bSJed Brown   TS          ts;                     /* timestepping context */
83c4762a1bSJed Brown   Mat         A;                      /* matrix data structure */
84c4762a1bSJed Brown   Vec         u;                      /* approximate solution vector */
85c4762a1bSJed Brown   PetscReal   time_total_max = 100.0; /* default max total time */
86c4762a1bSJed Brown   PetscInt    time_steps_max = 100;   /* default max timesteps */
87c4762a1bSJed Brown   PetscDraw   draw;                   /* drawing context */
88c4762a1bSJed Brown   PetscInt    steps, m;
89c4762a1bSJed Brown   PetscMPIInt size;
90c4762a1bSJed Brown   PetscReal   dt;
91c4762a1bSJed Brown   PetscBool   flg, flg_string;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown      Initialize program and set problem parameters
95c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96c4762a1bSJed Brown 
97327415f7SBarry Smith   PetscFunctionBeginUser;
98c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1003c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   m = 60;
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
105c4762a1bSJed Brown   flg_string = PETSC_FALSE;
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_string_viewer", &flg_string, NULL));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   appctx.m        = m;
109c4762a1bSJed Brown   appctx.h        = 1.0 / (m - 1.0);
110c4762a1bSJed Brown   appctx.norm_2   = 0.0;
111c4762a1bSJed Brown   appctx.norm_max = 0.0;
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving a linear TS problem on 1 processor\n"));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown      Create vector data structures
117c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /*
120c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
121c4762a1bSJed Brown   */
1229566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &u));
1239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 380, 400, 160, &appctx.viewer1));
1309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer1, 0, &draw));
1319566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
1329566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 0, 400, 160, &appctx.viewer2));
1339566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer2, 0, &draw));
1349566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown      Create timestepping solver context
138c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1419566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_LINEAR));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown      Set optional user-defined monitoring routine
145c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146c4762a1bSJed Brown 
14748a46eb9SPierre Jolivet   if (!flg_string) PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown 
151c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
152c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
1559566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m));
1569566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1579566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   flg = PETSC_FALSE;
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_ifunc", &flg, NULL));
161c4762a1bSJed Brown   if (!flg) {
162c4762a1bSJed Brown     appctx.A = NULL;
1639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-time_dependent_rhs", &flg, NULL));
164c4762a1bSJed Brown     if (flg) {
165c4762a1bSJed Brown       /*
166c4762a1bSJed Brown          For linear problems with a time-dependent f(u,t) in the equation
167dd8e379bSPierre Jolivet          u_t = f(u,t), the user provides the discretized right-hand side
168c4762a1bSJed Brown          as a time-dependent matrix.
169c4762a1bSJed Brown       */
1709566063dSJacob Faibussowitsch       PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1719566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(ts, A, A, RHSMatrixHeat, &appctx));
172c4762a1bSJed Brown     } else {
173c4762a1bSJed Brown       /*
174c4762a1bSJed Brown          For linear problems with a time-independent f(u) in the equation
175dd8e379bSPierre Jolivet          u_t = f(u), the user provides the discretized right-hand side
176c4762a1bSJed Brown          as a matrix only once, and then sets the special Jacobian evaluation
177c4762a1bSJed Brown          routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian.
178c4762a1bSJed Brown       */
1799566063dSJacob Faibussowitsch       PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx));
1809566063dSJacob Faibussowitsch       PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1819566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, &appctx));
182c4762a1bSJed Brown     }
183c4762a1bSJed Brown   } else {
184c4762a1bSJed Brown     Mat J;
185c4762a1bSJed Brown 
1869566063dSJacob Faibussowitsch     PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx));
1879566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &J));
1889566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunctionHeat, &appctx));
1899566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, J, J, IJacobianHeat, &appctx));
1909566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&J));
191c4762a1bSJed Brown 
1929566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
193c4762a1bSJed Brown     appctx.A      = A;
194c4762a1bSJed Brown     appctx.oshift = PETSC_MIN_REAL;
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197c4762a1bSJed Brown      Set solution vector and initial timestep
198c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   dt = appctx.h * appctx.h / 2.0;
2019566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown      Customize timestepping solver:
205c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
206c4762a1bSJed Brown        - Set timestepping duration info
207c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
208c4762a1bSJed Brown      For example,
209c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
210c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
211c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212c4762a1bSJed Brown 
2139566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
2149566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
2159566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2169566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219c4762a1bSJed Brown      Solve the problem
220c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /*
223c4762a1bSJed Brown      Evaluate initial conditions
224c4762a1bSJed Brown   */
2259566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /*
228c4762a1bSJed Brown      Run the timestepping solver
229c4762a1bSJed Brown   */
2309566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
2319566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234c4762a1bSJed Brown      View timestepping solver info
235c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236c4762a1bSJed Brown 
2379566063dSJacob 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)));
238c4762a1bSJed Brown   if (!flg_string) {
2399566063dSJacob Faibussowitsch     PetscCall(TSView(ts, PETSC_VIEWER_STDOUT_SELF));
240c4762a1bSJed Brown   } else {
241c4762a1bSJed Brown     PetscViewer stringviewer;
242c4762a1bSJed Brown     char        string[512];
243c4762a1bSJed Brown     const char *outstring;
244c4762a1bSJed Brown 
2459566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringOpen(PETSC_COMM_WORLD, string, sizeof(string), &stringviewer));
2469566063dSJacob Faibussowitsch     PetscCall(TSView(ts, stringviewer));
2479566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringGetStringRead(stringviewer, &outstring, NULL));
2483c633725SBarry Smith     PetscCheck((char *)outstring == (char *)string, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "String returned from viewer does not equal original string");
2499566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output from string viewer:%s\n", outstring));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&stringviewer));
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
254c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
255c4762a1bSJed Brown      are no longer needed.
256c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
257c4762a1bSJed Brown 
2589566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2619566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer1));
2629566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer2));
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
2649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.A));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /*
267c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
268c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
269c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
270c4762a1bSJed Brown          options are chosen (e.g., -log_view).
271c4762a1bSJed Brown   */
2729566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
273b122ec5aSJacob Faibussowitsch   return 0;
274c4762a1bSJed Brown }
275c4762a1bSJed Brown /* --------------------------------------------------------------------- */
276c4762a1bSJed Brown /*
277c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
278c4762a1bSJed Brown 
279c4762a1bSJed Brown    Input Parameter:
280c4762a1bSJed Brown    u - uninitialized solution vector (global)
281c4762a1bSJed Brown    appctx - user-defined application context
282c4762a1bSJed Brown 
283c4762a1bSJed Brown    Output Parameter:
284c4762a1bSJed Brown    u - vector with solution at initial time (global)
285c4762a1bSJed Brown */
InitialConditions(Vec u,AppCtx * appctx)286d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
287d71ae5a4SJacob Faibussowitsch {
288c4762a1bSJed Brown   PetscScalar *u_localptr, h = appctx->h;
289c4762a1bSJed Brown   PetscInt     i;
290c4762a1bSJed Brown 
2913ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
292c4762a1bSJed Brown   /*
293c4762a1bSJed Brown     Get a pointer to vector data.
294c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
295c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
296c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
297c4762a1bSJed Brown       the array.
298c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
299c4762a1bSJed Brown       C version.  See the users manual for details.
300c4762a1bSJed Brown   */
3019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(u, &u_localptr));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   /*
304c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
305c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
306c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
307c4762a1bSJed Brown   */
308c4762a1bSJed Brown   for (i = 0; i < appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h);
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   /*
311c4762a1bSJed Brown      Restore vector
312c4762a1bSJed Brown   */
3139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(u, &u_localptr));
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   /*
316c4762a1bSJed Brown      Print debugging information if desired
317c4762a1bSJed Brown   */
318c4762a1bSJed Brown   if (appctx->debug) {
3199566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess vector\n"));
3209566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
321c4762a1bSJed Brown   }
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown /* --------------------------------------------------------------------- */
325c4762a1bSJed Brown /*
326c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
327c4762a1bSJed Brown 
328c4762a1bSJed Brown    Input Parameters:
329c4762a1bSJed Brown    t - current time
330c4762a1bSJed Brown    solution - vector in which exact solution will be computed
331c4762a1bSJed Brown    appctx - user-defined application context
332c4762a1bSJed Brown 
333c4762a1bSJed Brown    Output Parameter:
334c4762a1bSJed Brown    solution - vector with the newly computed exact solution
335c4762a1bSJed Brown */
ExactSolution(PetscReal t,Vec solution,AppCtx * appctx)336d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx)
337d71ae5a4SJacob Faibussowitsch {
338c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2, tc = t;
339c4762a1bSJed Brown   PetscInt     i;
340c4762a1bSJed Brown 
3413ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
342c4762a1bSJed Brown   /*
343c4762a1bSJed Brown      Get a pointer to vector data.
344c4762a1bSJed Brown   */
3459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(solution, &s_localptr));
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   /*
348c4762a1bSJed Brown      Simply write the solution directly into the array locations.
349c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
350c4762a1bSJed Brown   */
351c4762a1bSJed Brown   ex1 = PetscExpScalar(-36. * PETSC_PI * PETSC_PI * tc);
352c4762a1bSJed Brown   ex2 = PetscExpScalar(-4. * PETSC_PI * PETSC_PI * tc);
3539371c9d4SSatish Balay   sc1 = PETSC_PI * 6. * h;
3549371c9d4SSatish Balay   sc2 = PETSC_PI * 2. * h;
355c4762a1bSJed Brown   for (i = 0; i < appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1 * (PetscReal)i) * ex1 + 3. * PetscSinScalar(sc2 * (PetscReal)i) * ex2;
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /*
358c4762a1bSJed Brown      Restore vector
359c4762a1bSJed Brown   */
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(solution, &s_localptr));
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
362c4762a1bSJed Brown }
363c4762a1bSJed Brown /* --------------------------------------------------------------------- */
364c4762a1bSJed Brown /*
365c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
366c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
367c4762a1bSJed Brown    error in two different norms.
368c4762a1bSJed Brown 
369c4762a1bSJed Brown    This example also demonstrates changing the timestep via TSSetTimeStep().
370c4762a1bSJed Brown 
371c4762a1bSJed Brown    Input Parameters:
372c4762a1bSJed Brown    ts     - the timestep context
373c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
374c4762a1bSJed Brown              initial condition)
375c4762a1bSJed Brown    time   - the current time
376c4762a1bSJed Brown    u      - the solution at this timestep
377c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
378c4762a1bSJed Brown             In this case we use the application context which contains
379c4762a1bSJed Brown             information about the problem size, workspace and the exact
380c4762a1bSJed Brown             solution.
381c4762a1bSJed Brown */
Monitor(TS ts,PetscInt step,PetscReal time,Vec u,PetscCtx ctx)382*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, PetscCtx ctx)
383d71ae5a4SJacob Faibussowitsch {
384c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
385c4762a1bSJed Brown   PetscReal norm_2, norm_max, dt, dttol;
386c4762a1bSJed Brown 
3873ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
388c4762a1bSJed Brown   /*
389c4762a1bSJed Brown      View a graph of the current iterate
390c4762a1bSJed Brown   */
3919566063dSJacob Faibussowitsch   PetscCall(VecView(u, appctx->viewer2));
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   /*
394c4762a1bSJed Brown      Compute the exact solution
395c4762a1bSJed Brown   */
3969566063dSJacob Faibussowitsch   PetscCall(ExactSolution(time, appctx->solution, appctx));
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   /*
399c4762a1bSJed Brown      Print debugging information if desired
400c4762a1bSJed Brown   */
401c4762a1bSJed Brown   if (appctx->debug) {
4029566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Computed solution vector\n"));
4039566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
4049566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Exact solution vector\n"));
4059566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
406c4762a1bSJed Brown   }
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   /*
409c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
410c4762a1bSJed Brown   */
4119566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
4129566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2));
413c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h) * norm_2;
4149566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max));
415c4762a1bSJed Brown 
4169566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
41763a3b9bcSJacob 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));
418c4762a1bSJed Brown 
419c4762a1bSJed Brown   appctx->norm_2 += norm_2;
420c4762a1bSJed Brown   appctx->norm_max += norm_max;
421c4762a1bSJed Brown 
422c4762a1bSJed Brown   dttol = .0001;
4239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-dttol", &dttol, NULL));
424c4762a1bSJed Brown   if (dt < dttol) {
425c4762a1bSJed Brown     dt *= .999;
4269566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
427c4762a1bSJed Brown   }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   /*
430c4762a1bSJed Brown      View a graph of the error
431c4762a1bSJed Brown   */
4329566063dSJacob Faibussowitsch   PetscCall(VecView(appctx->solution, appctx->viewer1));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /*
435c4762a1bSJed Brown      Print debugging information if desired
436c4762a1bSJed Brown   */
437c4762a1bSJed Brown   if (appctx->debug) {
4389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error vector\n"));
4399566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
440c4762a1bSJed Brown   }
4413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
442c4762a1bSJed Brown }
443c4762a1bSJed Brown /* --------------------------------------------------------------------- */
444c4762a1bSJed Brown /*
445c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
446c4762a1bSJed Brown    matrix for the heat equation.
447c4762a1bSJed Brown 
448c4762a1bSJed Brown    Input Parameters:
449c4762a1bSJed Brown    ts - the TS context
450c4762a1bSJed Brown    t - current time
451c4762a1bSJed Brown    global_in - global input vector
452c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
453c4762a1bSJed Brown 
454c4762a1bSJed Brown    Output Parameters:
455c4762a1bSJed Brown    AA - Jacobian matrix
4567addb90fSBarry Smith    BB - optionally different matrix used to construct the preconditioner
457c4762a1bSJed Brown 
458c4762a1bSJed Brown    Notes:
459c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
460c4762a1bSJed Brown    in Fortran as well as in C.
461c4762a1bSJed Brown */
RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,PetscCtx ctx)462*2a8381b2SBarry Smith PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec X, Mat AA, Mat BB, PetscCtx ctx)
463d71ae5a4SJacob Faibussowitsch {
464c4762a1bSJed Brown   Mat         A      = AA;            /* Jacobian matrix */
465c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
466c4762a1bSJed Brown   PetscInt    mstart = 0;
467c4762a1bSJed Brown   PetscInt    mend   = appctx->m;
468c4762a1bSJed Brown   PetscInt    i, idx[3];
469c4762a1bSJed Brown   PetscScalar v[3], stwo = -2. / (appctx->h * appctx->h), sone = -.5 * stwo;
470c4762a1bSJed Brown 
4713ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
472c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
473c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
474c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
475c4762a1bSJed Brown   /*
476c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
477c4762a1bSJed Brown   */
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   mstart = 0;
480c4762a1bSJed Brown   v[0]   = 1.0;
4819566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
482c4762a1bSJed Brown   mstart++;
483c4762a1bSJed Brown 
484c4762a1bSJed Brown   mend--;
485c4762a1bSJed Brown   v[0] = 1.0;
4869566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mend, 1, &mend, v, INSERT_VALUES));
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   /*
489c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
490c4762a1bSJed Brown      matrix one row at a time.
491c4762a1bSJed Brown   */
4929371c9d4SSatish Balay   v[0] = sone;
4939371c9d4SSatish Balay   v[1] = stwo;
4949371c9d4SSatish Balay   v[2] = sone;
495c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
4969371c9d4SSatish Balay     idx[0] = i - 1;
4979371c9d4SSatish Balay     idx[1] = i;
4989371c9d4SSatish Balay     idx[2] = i + 1;
4999566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
500c4762a1bSJed Brown   }
501c4762a1bSJed Brown 
502c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
503c4762a1bSJed Brown      Complete the matrix assembly process and set some options
504c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
505c4762a1bSJed Brown   /*
506c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
507c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
508c4762a1bSJed Brown      Computations can be done while messages are in transition
509c4762a1bSJed Brown      by placing code between these two statements.
510c4762a1bSJed Brown   */
5119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
513c4762a1bSJed Brown 
514c4762a1bSJed Brown   /*
515c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
516c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
517c4762a1bSJed Brown   */
5189566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
5193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
520c4762a1bSJed Brown }
521c4762a1bSJed Brown 
IFunctionHeat(TS ts,PetscReal t,Vec X,Vec Xdot,Vec r,PetscCtx ctx)522*2a8381b2SBarry Smith PetscErrorCode IFunctionHeat(TS ts, PetscReal t, Vec X, Vec Xdot, Vec r, PetscCtx ctx)
523d71ae5a4SJacob Faibussowitsch {
524c4762a1bSJed Brown   AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
525c4762a1bSJed Brown 
5263ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
5279566063dSJacob Faibussowitsch   PetscCall(MatMult(appctx->A, X, r));
5289566063dSJacob Faibussowitsch   PetscCall(VecAYPX(r, -1.0, Xdot));
5293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
530c4762a1bSJed Brown }
531c4762a1bSJed Brown 
IJacobianHeat(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal s,Mat A,Mat B,PetscCtx ctx)532*2a8381b2SBarry Smith PetscErrorCode IJacobianHeat(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal s, Mat A, Mat B, PetscCtx ctx)
533d71ae5a4SJacob Faibussowitsch {
534c4762a1bSJed Brown   AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
535c4762a1bSJed Brown 
5363ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
5373ba16761SJacob Faibussowitsch   if (appctx->oshift == s) PetscFunctionReturn(PETSC_SUCCESS);
5389566063dSJacob Faibussowitsch   PetscCall(MatCopy(appctx->A, A, SAME_NONZERO_PATTERN));
5399566063dSJacob Faibussowitsch   PetscCall(MatScale(A, -1));
5409566063dSJacob Faibussowitsch   PetscCall(MatShift(A, s));
5419566063dSJacob Faibussowitsch   PetscCall(MatCopy(A, B, SAME_NONZERO_PATTERN));
542c4762a1bSJed Brown   appctx->oshift = s;
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
544c4762a1bSJed Brown }
545c4762a1bSJed Brown 
546c4762a1bSJed Brown /*TEST
547c4762a1bSJed Brown 
548c4762a1bSJed Brown     test:
549188af4bfSBarry Smith       args: -nox -ts_type ssp -ts_time_step 0.0005
550c4762a1bSJed Brown 
551c4762a1bSJed Brown     test:
552c4762a1bSJed Brown       suffix: 2
553188af4bfSBarry Smith       args: -nox -ts_type ssp -ts_time_step 0.0005 -time_dependent_rhs 1
554c4762a1bSJed Brown 
555c4762a1bSJed Brown     test:
556c4762a1bSJed Brown       suffix: 3
557c4762a1bSJed Brown       args: -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason
558c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
559c4762a1bSJed Brown       requires: !single
560c4762a1bSJed Brown 
561c4762a1bSJed Brown     test:
562c4762a1bSJed Brown       suffix: 4
563c4762a1bSJed Brown       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason
564c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
565c4762a1bSJed Brown 
566c4762a1bSJed Brown     test:
567c4762a1bSJed Brown       suffix: 5
568c4762a1bSJed Brown       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs
569c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
570c4762a1bSJed Brown 
571c4762a1bSJed Brown     test:
572c4762a1bSJed Brown       requires: !single
573c4762a1bSJed Brown       suffix: pod_guess
574188af4bfSBarry Smith       args: -nox -ts_type beuler -use_ifunc -ts_time_step 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason
575c4762a1bSJed Brown 
576c4762a1bSJed Brown     test:
577c4762a1bSJed Brown       requires: !single
578c4762a1bSJed Brown       suffix: pod_guess_Ainner
579188af4bfSBarry Smith       args: -nox -ts_type beuler -use_ifunc -ts_time_step 0.0005 -ksp_guess_type pod -ksp_guess_pod_Ainner -pc_type none -ksp_converged_reason
580c4762a1bSJed Brown 
581c4762a1bSJed Brown     test:
582c4762a1bSJed Brown       requires: !single
583c4762a1bSJed Brown       suffix: fischer_guess
584188af4bfSBarry Smith       args: -nox -ts_type beuler -use_ifunc -ts_time_step 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason
585c4762a1bSJed Brown 
586c4762a1bSJed Brown     test:
587c4762a1bSJed Brown       requires: !single
588c4762a1bSJed Brown       suffix: fischer_guess_2
589188af4bfSBarry Smith       args: -nox -ts_type beuler -use_ifunc -ts_time_step 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 2,10 -pc_type none -ksp_converged_reason
590c4762a1bSJed Brown 
591c4762a1bSJed Brown     test:
592c4762a1bSJed Brown       requires: !single
593cbb17d71SDavid Wells       suffix: fischer_guess_3
594188af4bfSBarry Smith       args: -nox -ts_type beuler -use_ifunc -ts_time_step 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 3,10 -pc_type none -ksp_converged_reason
595cbb17d71SDavid Wells 
596cbb17d71SDavid Wells     test:
597cbb17d71SDavid Wells       requires: !single
598c4762a1bSJed Brown       suffix: stringview
599c4762a1bSJed Brown       args: -nox -ts_type rosw -test_string_viewer
600c4762a1bSJed Brown 
601c4762a1bSJed Brown     test:
602c4762a1bSJed Brown       requires: !single
603c4762a1bSJed Brown       suffix: stringview_euler
604c4762a1bSJed Brown       args: -nox -ts_type euler -test_string_viewer
605c4762a1bSJed Brown 
606c4762a1bSJed Brown TEST*/
607