xref: /petsc/src/ts/tutorials/ex2.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves a time-dependent nonlinear PDE. Uses implicit\n\
3c4762a1bSJed Brown timestepping.  Runtime options include:\n\
4c4762a1bSJed Brown   -M <xg>, where <xg> = number of grid points\n\
5c4762a1bSJed Brown   -debug : Activate debugging printouts\n\
6c4762a1bSJed Brown   -nox   : Deactivate x-window graphics\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /* ------------------------------------------------------------------------
9c4762a1bSJed Brown 
10c4762a1bSJed Brown    This program solves the PDE
11c4762a1bSJed Brown 
12c4762a1bSJed Brown                u * u_xx
13c4762a1bSJed Brown          u_t = ---------
14c4762a1bSJed Brown                2*(t+1)^2
15c4762a1bSJed Brown 
16c4762a1bSJed Brown     on the domain 0 <= x <= 1, with boundary conditions
17c4762a1bSJed Brown          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
18c4762a1bSJed Brown     and initial condition
19c4762a1bSJed Brown          u(0,x) = 1 + x*x.
20c4762a1bSJed Brown 
21c4762a1bSJed Brown     The exact solution is:
22c4762a1bSJed Brown          u(t,x) = (1 + x*x) * (1 + t)
23c4762a1bSJed Brown 
24c4762a1bSJed Brown     Note that since the solution is linear in time and quadratic in x,
25c4762a1bSJed Brown     the finite difference scheme actually computes the "exact" solution.
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     We use by default the backward Euler method.
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   ------------------------------------------------------------------------- */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown    Include "petscts.h" to use the PETSc timestepping routines. Note that
33c4762a1bSJed Brown    this file automatically includes "petscsys.h" and other lower-level
34c4762a1bSJed Brown    PETSc include files.
35c4762a1bSJed Brown 
36c4762a1bSJed Brown    Include the "petscdmda.h" to allow us to use the distributed array data
37c4762a1bSJed Brown    structures to manage the parallel grid.
38c4762a1bSJed Brown */
39c4762a1bSJed Brown #include <petscts.h>
40c4762a1bSJed Brown #include <petscdm.h>
41c4762a1bSJed Brown #include <petscdmda.h>
42c4762a1bSJed Brown #include <petscdraw.h>
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*
45c4762a1bSJed Brown    User-defined application context - contains data needed by the
46c4762a1bSJed Brown    application-provided callback routines.
47c4762a1bSJed Brown */
48c4762a1bSJed Brown typedef struct {
49c4762a1bSJed Brown   MPI_Comm  comm;      /* communicator */
50c4762a1bSJed Brown   DM        da;        /* distributed array data structure */
51c4762a1bSJed Brown   Vec       localwork; /* local ghosted work vector */
52c4762a1bSJed Brown   Vec       u_local;   /* local ghosted approximate solution vector */
53c4762a1bSJed Brown   Vec       solution;  /* global exact solution vector */
54c4762a1bSJed Brown   PetscInt  m;         /* total number of grid points */
55c4762a1bSJed Brown   PetscReal h;         /* mesh width: h = 1/(m-1) */
56c4762a1bSJed Brown   PetscBool debug;     /* flag (1 indicates activation of debugging printouts) */
57c4762a1bSJed Brown } AppCtx;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*
60c4762a1bSJed Brown    User-defined routines, provided below.
61c4762a1bSJed Brown */
62c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
63c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
64c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
65c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
66c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *);
67c4762a1bSJed Brown 
689371c9d4SSatish Balay int main(int argc, char **argv) {
69c4762a1bSJed Brown   AppCtx    appctx;               /* user-defined application context */
70c4762a1bSJed Brown   TS        ts;                   /* timestepping context */
71c4762a1bSJed Brown   Mat       A;                    /* Jacobian matrix data structure */
72c4762a1bSJed Brown   Vec       u;                    /* approximate solution vector */
73c4762a1bSJed Brown   PetscInt  time_steps_max = 100; /* default max timesteps */
74c4762a1bSJed Brown   PetscReal dt;
75c4762a1bSJed Brown   PetscReal time_total_max = 100.0; /* default max total time */
76c4762a1bSJed Brown   PetscBool mymonitor      = PETSC_FALSE;
77c4762a1bSJed Brown   PetscReal bounds[]       = {1.0, 3.3};
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80c4762a1bSJed Brown      Initialize program and set problem parameters
81c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82c4762a1bSJed Brown 
83327415f7SBarry Smith   PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
859566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, bounds));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
88c4762a1bSJed Brown   appctx.m    = 60;
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &appctx.m, NULL));
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-mymonitor", &mymonitor));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   appctx.h = 1.0 / (appctx.m - 1.0);
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Create vector data structures
98c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*
101c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
102c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
103c4762a1bSJed Brown      total grid values spread equally among all the processors.
104c4762a1bSJed Brown   */
1059566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, appctx.m, 1, 1, NULL, &appctx.da));
1069566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(appctx.da));
1079566063dSJacob Faibussowitsch   PetscCall(DMSetUp(appctx.da));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /*
110c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
111c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
112c4762a1bSJed Brown      have the same types.
113c4762a1bSJed Brown   */
1149566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(appctx.da, &u));
1159566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(appctx.da, &appctx.u_local));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /*
118c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
119c4762a1bSJed Brown      create global work vector for storing exact solution.
120c4762a1bSJed Brown   */
1219566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(appctx.u_local, &appctx.localwork));
1229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125c4762a1bSJed Brown      Create timestepping solver context; set callback routine for
126c4762a1bSJed Brown      right-hand-side function evaluation.
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1309566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1319566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &appctx));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134c4762a1bSJed Brown      Set optional user-defined monitoring routine
135c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136c4762a1bSJed Brown 
137*48a46eb9SPierre Jolivet   if (mymonitor) PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      For nonlinear problems, the user can provide a Jacobian evaluation
141c4762a1bSJed Brown      routine (or use a finite differencing approximation).
142c4762a1bSJed Brown 
143c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine.
144c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1479566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, appctx.m, appctx.m));
1489566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1499566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1509566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &appctx));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown      Set solution vector and initial timestep
154c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   dt = appctx.h / 2.0;
1579566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160c4762a1bSJed Brown      Customize timestepping solver:
161c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
162c4762a1bSJed Brown        - Set timestepping duration info
163c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
164c4762a1bSJed Brown      For example,
165c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
166c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
167c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
1709566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
1719566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
1729566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1739566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176c4762a1bSJed Brown      Solve the problem
177c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /*
180c4762a1bSJed Brown      Evaluate initial conditions
181c4762a1bSJed Brown   */
1829566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /*
185c4762a1bSJed Brown      Run the timestepping solver
186c4762a1bSJed Brown   */
1879566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
191c4762a1bSJed Brown      are no longer needed.
192c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1979566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&appctx.da));
1989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.localwork));
1999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
2009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.u_local));
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /*
203c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
204c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
205c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
206c4762a1bSJed Brown          options are chosen (e.g., -log_view).
207c4762a1bSJed Brown   */
2089566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
209b122ec5aSJacob Faibussowitsch   return 0;
210c4762a1bSJed Brown }
211c4762a1bSJed Brown /* --------------------------------------------------------------------- */
212c4762a1bSJed Brown /*
213c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
214c4762a1bSJed Brown 
215c4762a1bSJed Brown    Input Parameters:
216c4762a1bSJed Brown    u - uninitialized solution vector (global)
217c4762a1bSJed Brown    appctx - user-defined application context
218c4762a1bSJed Brown 
219c4762a1bSJed Brown    Output Parameter:
220c4762a1bSJed Brown    u - vector with solution at initial time (global)
221c4762a1bSJed Brown */
2229371c9d4SSatish Balay PetscErrorCode InitialConditions(Vec u, AppCtx *appctx) {
223c4762a1bSJed Brown   PetscScalar *u_localptr, h = appctx->h, x;
224c4762a1bSJed Brown   PetscInt     i, mybase, myend;
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /*
227c4762a1bSJed Brown      Determine starting point of each processor's range of
228c4762a1bSJed Brown      grid values.
229c4762a1bSJed Brown   */
2309566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(u, &mybase, &myend));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /*
233c4762a1bSJed Brown     Get a pointer to vector data.
234c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
235c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
236c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
237c4762a1bSJed Brown       the array.
238c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
239c4762a1bSJed Brown       C version.  See the users manual for details.
240c4762a1bSJed Brown   */
2419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &u_localptr));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /*
244c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
245c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
246c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
247c4762a1bSJed Brown   */
248c4762a1bSJed Brown   for (i = mybase; i < myend; i++) {
249c4762a1bSJed Brown     x                      = h * (PetscReal)i; /* current location in global grid */
250c4762a1bSJed Brown     u_localptr[i - mybase] = 1.0 + x * x;
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   /*
254c4762a1bSJed Brown      Restore vector
255c4762a1bSJed Brown   */
2569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &u_localptr));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   /*
259c4762a1bSJed Brown      Print debugging information if desired
260c4762a1bSJed Brown   */
261c4762a1bSJed Brown   if (appctx->debug) {
2629566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "initial guess vector\n"));
2639566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
264c4762a1bSJed Brown   }
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   return 0;
267c4762a1bSJed Brown }
268c4762a1bSJed Brown /* --------------------------------------------------------------------- */
269c4762a1bSJed Brown /*
270c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    Input Parameters:
273c4762a1bSJed Brown    t - current time
274c4762a1bSJed Brown    solution - vector in which exact solution will be computed
275c4762a1bSJed Brown    appctx - user-defined application context
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    Output Parameter:
278c4762a1bSJed Brown    solution - vector with the newly computed exact solution
279c4762a1bSJed Brown */
2809371c9d4SSatish Balay PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx) {
281c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, x;
282c4762a1bSJed Brown   PetscInt     i, mybase, myend;
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /*
285c4762a1bSJed Brown      Determine starting and ending points of each processor's
286c4762a1bSJed Brown      range of grid values
287c4762a1bSJed Brown   */
2889566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(solution, &mybase, &myend));
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   /*
291c4762a1bSJed Brown      Get a pointer to vector data.
292c4762a1bSJed Brown   */
2939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(solution, &s_localptr));
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /*
296c4762a1bSJed Brown      Simply write the solution directly into the array locations.
297c4762a1bSJed Brown      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
298c4762a1bSJed Brown   */
299c4762a1bSJed Brown   for (i = mybase; i < myend; i++) {
300c4762a1bSJed Brown     x                      = h * (PetscReal)i;
301c4762a1bSJed Brown     s_localptr[i - mybase] = (t + 1.0) * (1.0 + x * x);
302c4762a1bSJed Brown   }
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /*
305c4762a1bSJed Brown      Restore vector
306c4762a1bSJed Brown   */
3079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(solution, &s_localptr));
308c4762a1bSJed Brown   return 0;
309c4762a1bSJed Brown }
310c4762a1bSJed Brown /* --------------------------------------------------------------------- */
311c4762a1bSJed Brown /*
312c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
313c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
314c4762a1bSJed Brown    error in two different norms.
315c4762a1bSJed Brown 
316c4762a1bSJed Brown    Input Parameters:
317c4762a1bSJed Brown    ts     - the timestep context
318c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
319c4762a1bSJed Brown             initial condition)
320c4762a1bSJed Brown    time   - the current time
321c4762a1bSJed Brown    u      - the solution at this timestep
322c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
323c4762a1bSJed Brown             In this case we use the application context which contains
324c4762a1bSJed Brown             information about the problem size, workspace and the exact
325c4762a1bSJed Brown             solution.
326c4762a1bSJed Brown */
3279371c9d4SSatish Balay PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) {
328c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
329c4762a1bSJed Brown   PetscReal en2, en2s, enmax;
330c4762a1bSJed Brown   PetscDraw draw;
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   /*
333e1dfdf8eSBarry Smith      We use the default X Windows viewer
334c4762a1bSJed Brown              PETSC_VIEWER_DRAW_(appctx->comm)
335c4762a1bSJed Brown      that is associated with the current communicator. This saves
336c4762a1bSJed Brown      the effort of calling PetscViewerDrawOpen() to create the window.
337c4762a1bSJed Brown      Note that if we wished to plot several items in separate windows we
338c4762a1bSJed Brown      would create each viewer with PetscViewerDrawOpen() and store them in
339c4762a1bSJed Brown      the application context, appctx.
340c4762a1bSJed Brown 
341c4762a1bSJed Brown      PetscReal buffering makes graphics look better.
342c4762a1bSJed Brown   */
3439566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm), 0, &draw));
3449566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
3459566063dSJacob Faibussowitsch   PetscCall(VecView(u, PETSC_VIEWER_DRAW_(appctx->comm)));
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   /*
348c4762a1bSJed Brown      Compute the exact solution at this timestep
349c4762a1bSJed Brown   */
3509566063dSJacob Faibussowitsch   PetscCall(ExactSolution(time, appctx->solution, appctx));
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   /*
353c4762a1bSJed Brown      Print debugging information if desired
354c4762a1bSJed Brown   */
355c4762a1bSJed Brown   if (appctx->debug) {
3569566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Computed solution vector\n"));
3579566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
3589566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Exact solution vector\n"));
3599566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD));
360c4762a1bSJed Brown   }
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   /*
363c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
364c4762a1bSJed Brown   */
3659566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
3669566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &en2));
367c4762a1bSJed Brown   en2s = PetscSqrtReal(appctx->h) * en2; /* scale the 2-norm by the grid spacing */
3689566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &enmax));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   /*
371c4762a1bSJed Brown      PetscPrintf() causes only the first processor in this
372c4762a1bSJed Brown      communicator to print the timestep information.
373c4762a1bSJed Brown   */
37463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(appctx->comm, "Timestep %" PetscInt_FMT ": time = %g 2-norm error = %g  max norm error = %g\n", step, (double)time, (double)en2s, (double)enmax));
375c4762a1bSJed Brown 
376c4762a1bSJed Brown   /*
377c4762a1bSJed Brown      Print debugging information if desired
378c4762a1bSJed Brown   */
379c4762a1bSJed Brown   if (appctx->debug) {
3809566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Error vector\n"));
3819566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD));
382c4762a1bSJed Brown   }
383c4762a1bSJed Brown   return 0;
384c4762a1bSJed Brown }
385c4762a1bSJed Brown /* --------------------------------------------------------------------- */
386c4762a1bSJed Brown /*
387c4762a1bSJed Brown    RHSFunction - User-provided routine that evalues the right-hand-side
388c4762a1bSJed Brown    function of the ODE.  This routine is set in the main program by
389c4762a1bSJed Brown    calling TSSetRHSFunction().  We compute:
390c4762a1bSJed Brown           global_out = F(global_in)
391c4762a1bSJed Brown 
392c4762a1bSJed Brown    Input Parameters:
393c4762a1bSJed Brown    ts         - timesteping context
394c4762a1bSJed Brown    t          - current time
395c4762a1bSJed Brown    global_in  - vector containing the current iterate
396c4762a1bSJed Brown    ctx        - (optional) user-provided context for function evaluation.
397c4762a1bSJed Brown                 In this case we use the appctx defined above.
398c4762a1bSJed Brown 
399c4762a1bSJed Brown    Output Parameter:
400c4762a1bSJed Brown    global_out - vector containing the newly evaluated function
401c4762a1bSJed Brown */
4029371c9d4SSatish Balay PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec global_in, Vec global_out, void *ctx) {
403c4762a1bSJed Brown   AppCtx            *appctx    = (AppCtx *)ctx;     /* user-defined application context */
404c4762a1bSJed Brown   DM                 da        = appctx->da;        /* distributed array */
405c4762a1bSJed Brown   Vec                local_in  = appctx->u_local;   /* local ghosted input vector */
406c4762a1bSJed Brown   Vec                localwork = appctx->localwork; /* local ghosted work vector */
407c4762a1bSJed Brown   PetscInt           i, localsize;
408c4762a1bSJed Brown   PetscMPIInt        rank, size;
409c4762a1bSJed Brown   PetscScalar       *copyptr, sc;
410c4762a1bSJed Brown   const PetscScalar *localptr;
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
413c4762a1bSJed Brown      Get ready for local function computations
414c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
415c4762a1bSJed Brown   /*
416c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
417c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
418c4762a1bSJed Brown      By placing code between these two statements, computations can be
419c4762a1bSJed Brown      done while messages are in transition.
420c4762a1bSJed Brown   */
4219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, global_in, INSERT_VALUES, local_in));
4229566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, global_in, INSERT_VALUES, local_in));
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /*
425c4762a1bSJed Brown       Access directly the values in our local INPUT work array
426c4762a1bSJed Brown   */
4279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in, &localptr));
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   /*
430c4762a1bSJed Brown       Access directly the values in our local OUTPUT work array
431c4762a1bSJed Brown   */
4329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localwork, &copyptr));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   sc = 1.0 / (appctx->h * appctx->h * 2.0 * (1.0 + t) * (1.0 + t));
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /*
437c4762a1bSJed Brown       Evaluate our function on the nodes owned by this processor
438c4762a1bSJed Brown   */
4399566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(local_in, &localsize));
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
442c4762a1bSJed Brown      Compute entries for the locally owned part
443c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /*
446c4762a1bSJed Brown      Handle boundary conditions: This is done by using the boundary condition
447c4762a1bSJed Brown         u(t,boundary) = g(t,boundary)
448c4762a1bSJed Brown      for some function g. Now take the derivative with respect to t to obtain
449c4762a1bSJed Brown         u_{t}(t,boundary) = g_{t}(t,boundary)
450c4762a1bSJed Brown 
451c4762a1bSJed Brown      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
452c4762a1bSJed Brown              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
453c4762a1bSJed Brown   */
4549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(appctx->comm, &rank));
4559566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(appctx->comm, &size));
456dd400576SPatrick Sanan   if (rank == 0) copyptr[0] = 1.0;
457c4762a1bSJed Brown   if (rank == size - 1) copyptr[localsize - 1] = 2.0;
458c4762a1bSJed Brown 
459c4762a1bSJed Brown   /*
460c4762a1bSJed Brown      Handle the interior nodes where the PDE is replace by finite
461c4762a1bSJed Brown      difference operators.
462c4762a1bSJed Brown   */
463c4762a1bSJed Brown   for (i = 1; i < localsize - 1; i++) copyptr[i] = localptr[i] * sc * (localptr[i + 1] + localptr[i - 1] - 2.0 * localptr[i]);
464c4762a1bSJed Brown 
465c4762a1bSJed Brown   /*
466c4762a1bSJed Brown      Restore vectors
467c4762a1bSJed Brown   */
4689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in, &localptr));
4699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localwork, &copyptr));
470c4762a1bSJed Brown 
471c4762a1bSJed Brown   /*
472c4762a1bSJed Brown      Insert values from the local OUTPUT vector into the global
473c4762a1bSJed Brown      output vector
474c4762a1bSJed Brown   */
4759566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, localwork, INSERT_VALUES, global_out));
4769566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, localwork, INSERT_VALUES, global_out));
477c4762a1bSJed Brown 
478c4762a1bSJed Brown   /* Print debugging information if desired */
479c4762a1bSJed Brown   if (appctx->debug) {
4809566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "RHS function vector\n"));
4819566063dSJacob Faibussowitsch     PetscCall(VecView(global_out, PETSC_VIEWER_STDOUT_WORLD));
482c4762a1bSJed Brown   }
483c4762a1bSJed Brown 
484c4762a1bSJed Brown   return 0;
485c4762a1bSJed Brown }
486c4762a1bSJed Brown /* --------------------------------------------------------------------- */
487c4762a1bSJed Brown /*
488c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
489c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
490c4762a1bSJed Brown 
491c4762a1bSJed Brown    Input Parameters:
492c4762a1bSJed Brown    ts - the TS context
493c4762a1bSJed Brown    t - current time
494c4762a1bSJed Brown    global_in - global input vector
495c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
496c4762a1bSJed Brown 
497c4762a1bSJed Brown    Output Parameters:
498c4762a1bSJed Brown    AA - Jacobian matrix
499c4762a1bSJed Brown    BB - optionally different preconditioning matrix
500c4762a1bSJed Brown    str - flag indicating matrix structure
501c4762a1bSJed Brown 
502c4762a1bSJed Brown   Notes:
503c4762a1bSJed Brown   RHSJacobian computes entries for the locally owned part of the Jacobian.
504c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
505c4762a1bSJed Brown      contiguous chunks of rows across the processors.
506c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
507c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
508c4762a1bSJed Brown      appropriate processor during matrix assembly).
509c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
510c4762a1bSJed Brown      using MatSetValues().
511c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
512c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
513c4762a1bSJed Brown      in Fortran as well as in C.
514c4762a1bSJed Brown */
5159371c9d4SSatish Balay PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec global_in, Mat AA, Mat BB, void *ctx) {
516c4762a1bSJed Brown   AppCtx            *appctx   = (AppCtx *)ctx;   /* user-defined application context */
517c4762a1bSJed Brown   Vec                local_in = appctx->u_local; /* local ghosted input vector */
518c4762a1bSJed Brown   DM                 da       = appctx->da;      /* distributed array */
519c4762a1bSJed Brown   PetscScalar        v[3], sc;
520c4762a1bSJed Brown   const PetscScalar *localptr;
521c4762a1bSJed Brown   PetscInt           i, mstart, mend, mstarts, mends, idx[3], is;
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
524c4762a1bSJed Brown      Get ready for local Jacobian computations
525c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
526c4762a1bSJed Brown   /*
527c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
528c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
529c4762a1bSJed Brown      By placing code between these two statements, computations can be
530c4762a1bSJed Brown      done while messages are in transition.
531c4762a1bSJed Brown   */
5329566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, global_in, INSERT_VALUES, local_in));
5339566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, global_in, INSERT_VALUES, local_in));
534c4762a1bSJed Brown 
535c4762a1bSJed Brown   /*
536c4762a1bSJed Brown      Get pointer to vector data
537c4762a1bSJed Brown   */
5389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in, &localptr));
539c4762a1bSJed Brown 
540c4762a1bSJed Brown   /*
541c4762a1bSJed Brown      Get starting and ending locally owned rows of the matrix
542c4762a1bSJed Brown   */
5439566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(BB, &mstarts, &mends));
5449371c9d4SSatish Balay   mstart = mstarts;
5459371c9d4SSatish Balay   mend   = mends;
546c4762a1bSJed Brown 
547c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
548c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
549c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
550c4762a1bSJed Brown         contiguous chunks of rows across the processors.
551c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
552c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
553c4762a1bSJed Brown         appropriate processor during matrix assembly).
554c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
555c4762a1bSJed Brown       - We can set matrix entries either using either
556c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues().
557c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
558c4762a1bSJed Brown 
559c4762a1bSJed Brown   /*
560c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
561c4762a1bSJed Brown   */
562c4762a1bSJed Brown   if (mstart == 0) {
563c4762a1bSJed Brown     v[0] = 0.0;
5649566063dSJacob Faibussowitsch     PetscCall(MatSetValues(BB, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
565c4762a1bSJed Brown     mstart++;
566c4762a1bSJed Brown   }
567c4762a1bSJed Brown   if (mend == appctx->m) {
568c4762a1bSJed Brown     mend--;
569c4762a1bSJed Brown     v[0] = 0.0;
5709566063dSJacob Faibussowitsch     PetscCall(MatSetValues(BB, 1, &mend, 1, &mend, v, INSERT_VALUES));
571c4762a1bSJed Brown   }
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   /*
574c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
575c4762a1bSJed Brown      matrix one row at a time.
576c4762a1bSJed Brown   */
577c4762a1bSJed Brown   sc = 1.0 / (appctx->h * appctx->h * 2.0 * (1.0 + t) * (1.0 + t));
578c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
5799371c9d4SSatish Balay     idx[0] = i - 1;
5809371c9d4SSatish Balay     idx[1] = i;
5819371c9d4SSatish Balay     idx[2] = i + 1;
582c4762a1bSJed Brown     is     = i - mstart + 1;
583c4762a1bSJed Brown     v[0]   = sc * localptr[is];
584c4762a1bSJed Brown     v[1]   = sc * (localptr[is + 1] + localptr[is - 1] - 4.0 * localptr[is]);
585c4762a1bSJed Brown     v[2]   = sc * localptr[is];
5869566063dSJacob Faibussowitsch     PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES));
587c4762a1bSJed Brown   }
588c4762a1bSJed Brown 
589c4762a1bSJed Brown   /*
590c4762a1bSJed Brown      Restore vector
591c4762a1bSJed Brown   */
5929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in, &localptr));
593c4762a1bSJed Brown 
594c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
595c4762a1bSJed Brown      Complete the matrix assembly process and set some options
596c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
597c4762a1bSJed Brown   /*
598c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
599c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
600c4762a1bSJed Brown      Computations can be done while messages are in transition
601c4762a1bSJed Brown      by placing code between these two statements.
602c4762a1bSJed Brown   */
6039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY));
6049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY));
605c4762a1bSJed Brown   if (BB != AA) {
6069566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(AA, MAT_FINAL_ASSEMBLY));
6079566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(AA, MAT_FINAL_ASSEMBLY));
608c4762a1bSJed Brown   }
609c4762a1bSJed Brown 
610c4762a1bSJed Brown   /*
611c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
612c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
613c4762a1bSJed Brown   */
6149566063dSJacob Faibussowitsch   PetscCall(MatSetOption(BB, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
615c4762a1bSJed Brown 
616c4762a1bSJed Brown   return 0;
617c4762a1bSJed Brown }
618c4762a1bSJed Brown 
619c4762a1bSJed Brown /*TEST
620c4762a1bSJed Brown 
621c4762a1bSJed Brown     test:
622c4762a1bSJed Brown       args: -nox -ts_dt 10 -mymonitor
623c4762a1bSJed Brown       nsize: 2
624c4762a1bSJed Brown       requires: !single
625c4762a1bSJed Brown 
626c4762a1bSJed Brown     test:
627c4762a1bSJed Brown       suffix: tut_1
628c4762a1bSJed Brown       nsize: 1
629c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor
630c4762a1bSJed Brown 
631c4762a1bSJed Brown     test:
632c4762a1bSJed Brown       suffix: tut_2
633c4762a1bSJed Brown       nsize: 4
634c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -snes_monitor -ksp_monitor
635c4762a1bSJed Brown 
636c4762a1bSJed Brown     test:
637c4762a1bSJed Brown       suffix: tut_3
638c4762a1bSJed Brown       nsize: 4
6392e16c0ceSBarry Smith       args: -ts_max_steps 10 -ts_monitor -M 128
640c4762a1bSJed Brown 
641c4762a1bSJed Brown TEST*/
642