xref: /petsc/src/ts/tutorials/ex6.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   -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 
76d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
77d71ae5a4SJacob Faibussowitsch {
78c4762a1bSJed Brown   AppCtx      appctx;                 /* user-defined application context */
79c4762a1bSJed Brown   TS          ts;                     /* timestepping context */
80c4762a1bSJed Brown   Mat         A;                      /* matrix data structure */
81c4762a1bSJed Brown   Vec         u;                      /* approximate solution vector */
82c4762a1bSJed Brown   PetscReal   time_total_max = 100.0; /* default max total time */
83c4762a1bSJed Brown   PetscInt    time_steps_max = 100;   /* default max timesteps */
84c4762a1bSJed Brown   PetscDraw   draw;                   /* drawing context */
85c4762a1bSJed Brown   PetscInt    steps, m;
86c4762a1bSJed Brown   PetscMPIInt size;
87c4762a1bSJed Brown   PetscReal   dt;
88c4762a1bSJed Brown   PetscReal   ftime;
89c4762a1bSJed Brown   PetscBool   flg;
90c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91c4762a1bSJed Brown      Initialize program and set problem parameters
92c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93c4762a1bSJed Brown 
94327415f7SBarry Smith   PetscFunctionBeginUser;
959566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
96c4762a1bSJed Brown   MPI_Comm_size(PETSC_COMM_WORLD, &size);
973c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   m = 60;
1009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   appctx.m        = m;
104c4762a1bSJed Brown   appctx.h        = 1.0 / (m - 1.0);
105c4762a1bSJed Brown   appctx.norm_2   = 0.0;
106c4762a1bSJed Brown   appctx.norm_max = 0.0;
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving a linear TS problem on 1 processor\n"));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4762a1bSJed Brown      Create vector data structures
112c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /*
115c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
116c4762a1bSJed Brown   */
1179566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &u));
1189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
122c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 380, 400, 160, &appctx.viewer1));
1259566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer1, 0, &draw));
1269566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
1279566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 0, 400, 160, &appctx.viewer2));
1289566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer2, 0, &draw));
1299566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132c4762a1bSJed Brown      Create timestepping solver context
133c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134c4762a1bSJed Brown 
1359566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1369566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_LINEAR));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139c4762a1bSJed Brown      Set optional user-defined monitoring routine
140c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown 
146c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
147c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148c4762a1bSJed Brown 
1499566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
1509566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m));
1519566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1529566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-time_dependent_rhs", &flg));
155c4762a1bSJed Brown   if (flg) {
156c4762a1bSJed Brown     /*
157c4762a1bSJed Brown        For linear problems with a time-dependent f(u,t) in the equation
158c4762a1bSJed Brown        u_t = f(u,t), the user provides the discretized right-hand-side
159c4762a1bSJed Brown        as a time-dependent matrix.
160c4762a1bSJed Brown     */
1619566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1629566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, A, A, RHSMatrixHeat, &appctx));
163c4762a1bSJed Brown   } else {
164c4762a1bSJed Brown     /*
165c4762a1bSJed Brown        For linear problems with a time-independent f(u) in the equation
166c4762a1bSJed Brown        u_t = f(u), the user provides the discretized right-hand-side
167c4762a1bSJed Brown        as a matrix only once, and then sets a null matrix evaluation
168c4762a1bSJed Brown        routine.
169c4762a1bSJed Brown     */
1709566063dSJacob Faibussowitsch     PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx));
1719566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1729566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, &appctx));
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176c4762a1bSJed Brown      Set solution vector and initial timestep
177c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   dt = appctx.h * appctx.h / 2.0;
1809566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
1819566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, u));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184c4762a1bSJed Brown      Customize timestepping solver:
185c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
186c4762a1bSJed Brown        - Set timestepping duration info
187c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
188c4762a1bSJed Brown      For example,
189c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
190c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
191c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192c4762a1bSJed Brown 
1939566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
1949566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
1959566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1969566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199c4762a1bSJed Brown      Solve the problem
200c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /*
203c4762a1bSJed Brown      Evaluate initial conditions
204c4762a1bSJed Brown   */
2059566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /*
208c4762a1bSJed Brown      Run the timestepping solver
209c4762a1bSJed Brown   */
2109566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
2119566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2129566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215c4762a1bSJed Brown      View timestepping solver info
216c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217c4762a1bSJed Brown 
2189566063dSJacob 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)));
2199566063dSJacob Faibussowitsch   PetscCall(TSView(ts, PETSC_VIEWER_STDOUT_SELF));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
223c4762a1bSJed Brown      are no longer needed.
224c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2299566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer1));
2309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer2));
2319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /*
234c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
235c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
236c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
237c4762a1bSJed Brown          options are chosen (e.g., -log_view).
238c4762a1bSJed Brown   */
2399566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
240b122ec5aSJacob Faibussowitsch   return 0;
241c4762a1bSJed Brown }
242c4762a1bSJed Brown /* --------------------------------------------------------------------- */
243c4762a1bSJed Brown /*
244c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
245c4762a1bSJed Brown 
246c4762a1bSJed Brown    Input Parameter:
247c4762a1bSJed Brown    u - uninitialized solution vector (global)
248c4762a1bSJed Brown    appctx - user-defined application context
249c4762a1bSJed Brown 
250c4762a1bSJed Brown    Output Parameter:
251c4762a1bSJed Brown    u - vector with solution at initial time (global)
252c4762a1bSJed Brown */
253d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
254d71ae5a4SJacob Faibussowitsch {
255c4762a1bSJed Brown   PetscScalar *u_localptr;
256c4762a1bSJed Brown   PetscInt     i;
257c4762a1bSJed Brown 
258*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
259c4762a1bSJed Brown   /*
260c4762a1bSJed Brown     Get a pointer to vector data.
261c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
262c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
263c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
264c4762a1bSJed Brown       the array.
265c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
266c4762a1bSJed Brown       C version.  See the users manual for details.
267c4762a1bSJed Brown   */
2689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &u_localptr));
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   /*
271c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
272c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
273c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
274c4762a1bSJed Brown   */
275c4762a1bSJed 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);
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /*
278c4762a1bSJed Brown      Restore vector
279c4762a1bSJed Brown   */
2809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &u_localptr));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /*
283c4762a1bSJed Brown      Print debugging information if desired
284c4762a1bSJed Brown   */
2851baa6e33SBarry Smith   if (appctx->debug) PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
286c4762a1bSJed Brown 
287*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
288c4762a1bSJed Brown }
289c4762a1bSJed Brown /* --------------------------------------------------------------------- */
290c4762a1bSJed Brown /*
291c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    Input Parameters:
294c4762a1bSJed Brown    t - current time
295c4762a1bSJed Brown    solution - vector in which exact solution will be computed
296c4762a1bSJed Brown    appctx - user-defined application context
297c4762a1bSJed Brown 
298c4762a1bSJed Brown    Output Parameter:
299c4762a1bSJed Brown    solution - vector with the newly computed exact solution
300c4762a1bSJed Brown */
301d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx)
302d71ae5a4SJacob Faibussowitsch {
303c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
304c4762a1bSJed Brown   PetscInt     i;
305c4762a1bSJed Brown 
306*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
307c4762a1bSJed Brown   /*
308c4762a1bSJed Brown      Get a pointer to vector data.
309c4762a1bSJed Brown   */
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(solution, &s_localptr));
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /*
313c4762a1bSJed Brown      Simply write the solution directly into the array locations.
314c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
315c4762a1bSJed Brown   */
3169371c9d4SSatish Balay   ex1 = PetscExpReal(-36. * PETSC_PI * PETSC_PI * t);
3179371c9d4SSatish Balay   ex2 = PetscExpReal(-4. * PETSC_PI * PETSC_PI * t);
3189371c9d4SSatish Balay   sc1 = PETSC_PI * 6. * h;
3199371c9d4SSatish Balay   sc2 = PETSC_PI * 2. * h;
320c4762a1bSJed Brown   for (i = 0; i < appctx->m; i++) s_localptr[i] = PetscSinReal(PetscRealPart(sc1) * (PetscReal)i) * ex1 + 3. * PetscSinReal(PetscRealPart(sc2) * (PetscReal)i) * ex2;
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   /*
323c4762a1bSJed Brown      Restore vector
324c4762a1bSJed Brown   */
3259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(solution, &s_localptr));
326*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
327c4762a1bSJed Brown }
328c4762a1bSJed Brown /* --------------------------------------------------------------------- */
329c4762a1bSJed Brown /*
330c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
331c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
332c4762a1bSJed Brown    error in two different norms.
333c4762a1bSJed Brown 
334c4762a1bSJed Brown    This example also demonstrates changing the timestep via TSSetTimeStep().
335c4762a1bSJed Brown 
336c4762a1bSJed Brown    Input Parameters:
337c4762a1bSJed Brown    ts     - the timestep context
338c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
339c4762a1bSJed Brown              initial condition)
340c4762a1bSJed Brown    crtime  - the current time
341c4762a1bSJed Brown    u      - the solution at this timestep
342c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
343c4762a1bSJed Brown             In this case we use the application context which contains
344c4762a1bSJed Brown             information about the problem size, workspace and the exact
345c4762a1bSJed Brown             solution.
346c4762a1bSJed Brown */
347d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
348d71ae5a4SJacob Faibussowitsch {
349c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
350c4762a1bSJed Brown   PetscReal norm_2, norm_max, dt, dttol;
351c4762a1bSJed Brown   PetscBool flg;
352c4762a1bSJed Brown 
353*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
354c4762a1bSJed Brown   /*
355c4762a1bSJed Brown      View a graph of the current iterate
356c4762a1bSJed Brown   */
3579566063dSJacob Faibussowitsch   PetscCall(VecView(u, appctx->viewer2));
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   /*
360c4762a1bSJed Brown      Compute the exact solution
361c4762a1bSJed Brown   */
3629566063dSJacob Faibussowitsch   PetscCall(ExactSolution(crtime, appctx->solution, appctx));
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /*
365c4762a1bSJed Brown      Print debugging information if desired
366c4762a1bSJed Brown   */
367c4762a1bSJed Brown   if (appctx->debug) {
3689566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Computed solution vector\n"));
3699566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
3709566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Exact solution vector\n"));
3719566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
372c4762a1bSJed Brown   }
373c4762a1bSJed Brown 
374c4762a1bSJed Brown   /*
375c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
376c4762a1bSJed Brown   */
3779566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
3789566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2));
379c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h) * norm_2;
3809566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max));
381c4762a1bSJed Brown 
3829566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
38348a46eb9SPierre 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));
384c4762a1bSJed Brown   appctx->norm_2 += norm_2;
385c4762a1bSJed Brown   appctx->norm_max += norm_max;
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   dttol = .0001;
3889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-dttol", &dttol, &flg));
389c4762a1bSJed Brown   if (dt < dttol) {
390c4762a1bSJed Brown     dt *= .999;
3919566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
392c4762a1bSJed Brown   }
393c4762a1bSJed Brown 
394c4762a1bSJed Brown   /*
395c4762a1bSJed Brown      View a graph of the error
396c4762a1bSJed Brown   */
3979566063dSJacob Faibussowitsch   PetscCall(VecView(appctx->solution, appctx->viewer1));
398c4762a1bSJed Brown 
399c4762a1bSJed Brown   /*
400c4762a1bSJed Brown      Print debugging information if desired
401c4762a1bSJed Brown   */
402c4762a1bSJed Brown   if (appctx->debug) {
4039566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error vector\n"));
4049566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
405c4762a1bSJed Brown   }
406c4762a1bSJed Brown 
407*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408c4762a1bSJed Brown }
409c4762a1bSJed Brown /* --------------------------------------------------------------------- */
410c4762a1bSJed Brown /*
411c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
412c4762a1bSJed Brown    matrix for the heat equation.
413c4762a1bSJed Brown 
414c4762a1bSJed Brown    Input Parameters:
415c4762a1bSJed Brown    ts - the TS context
416c4762a1bSJed Brown    t - current time
417c4762a1bSJed Brown    global_in - global input vector
418c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
419c4762a1bSJed Brown 
420c4762a1bSJed Brown    Output Parameters:
421c4762a1bSJed Brown    AA - Jacobian matrix
422c4762a1bSJed Brown    BB - optionally different preconditioning matrix
423c4762a1bSJed Brown    str - flag indicating matrix structure
424c4762a1bSJed Brown 
425c4762a1bSJed Brown    Notes:
426c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
427c4762a1bSJed Brown    in Fortran as well as in C.
428c4762a1bSJed Brown */
429d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec X, Mat AA, Mat BB, void *ctx)
430d71ae5a4SJacob Faibussowitsch {
431c4762a1bSJed Brown   Mat         A      = AA;            /* Jacobian matrix */
432c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
433c4762a1bSJed Brown   PetscInt    mstart = 0;
434c4762a1bSJed Brown   PetscInt    mend   = appctx->m;
435c4762a1bSJed Brown   PetscInt    i, idx[3];
436c4762a1bSJed Brown   PetscScalar v[3], stwo = -2. / (appctx->h * appctx->h), sone = -.5 * stwo;
437c4762a1bSJed Brown 
438*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
439c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
440c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
441c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
442c4762a1bSJed Brown   /*
443c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
444c4762a1bSJed Brown   */
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   mstart = 0;
447c4762a1bSJed Brown   v[0]   = 1.0;
4489566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
449c4762a1bSJed Brown   mstart++;
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   mend--;
452c4762a1bSJed Brown   v[0] = 1.0;
4539566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &mend, 1, &mend, v, INSERT_VALUES));
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   /*
456c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
457c4762a1bSJed Brown      matrix one row at a time.
458c4762a1bSJed Brown   */
4599371c9d4SSatish Balay   v[0] = sone;
4609371c9d4SSatish Balay   v[1] = stwo;
4619371c9d4SSatish Balay   v[2] = sone;
462c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
4639371c9d4SSatish Balay     idx[0] = i - 1;
4649371c9d4SSatish Balay     idx[1] = i;
4659371c9d4SSatish Balay     idx[2] = i + 1;
4669566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
467c4762a1bSJed Brown   }
468c4762a1bSJed Brown 
469c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
470c4762a1bSJed Brown      Complete the matrix assembly process and set some options
471c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
472c4762a1bSJed Brown   /*
473c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
474c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
475c4762a1bSJed Brown      Computations can be done while messages are in transition
476c4762a1bSJed Brown      by placing code between these two statements.
477c4762a1bSJed Brown   */
4789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
480c4762a1bSJed Brown 
481c4762a1bSJed Brown   /*
482c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
483c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
484c4762a1bSJed Brown   */
4859566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
486c4762a1bSJed Brown 
487*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
488c4762a1bSJed Brown }
489c4762a1bSJed Brown /* --------------------------------------------------------------------- */
490c4762a1bSJed Brown /*
491c4762a1bSJed Brown    Input Parameters:
492c4762a1bSJed Brown    ts - the TS context
493c4762a1bSJed Brown    t - current time
494c4762a1bSJed Brown    f - function
495c4762a1bSJed Brown    ctx - optional user-defined context, as set by TSetBCFunction()
496c4762a1bSJed Brown  */
497d71ae5a4SJacob Faibussowitsch PetscErrorCode MyBCRoutine(TS ts, PetscReal t, Vec f, void *ctx)
498d71ae5a4SJacob Faibussowitsch {
499c4762a1bSJed Brown   AppCtx      *appctx = (AppCtx *)ctx; /* user-defined application context */
500c4762a1bSJed Brown   PetscInt     m      = appctx->m;
501c4762a1bSJed Brown   PetscScalar *fa;
502c4762a1bSJed Brown 
503*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
5049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &fa));
505c4762a1bSJed Brown   fa[0]     = 0.0;
506c4762a1bSJed Brown   fa[m - 1] = 1.0;
5079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &fa));
5089566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "t=%g\n", (double)t));
509c4762a1bSJed Brown 
510*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
511c4762a1bSJed Brown }
512c4762a1bSJed Brown 
513c4762a1bSJed Brown /*TEST
514c4762a1bSJed Brown 
515c4762a1bSJed Brown     test:
516c4762a1bSJed Brown       args: -nox -ts_max_steps 4
517c4762a1bSJed Brown 
518c4762a1bSJed Brown TEST*/
519