xref: /petsc/src/ts/tutorials/ex4.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   -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 one-dimensional heat equation (also called the
11c4762a1bSJed Brown    diffusion equation),
12c4762a1bSJed Brown        u_t = u_xx,
13c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
14c4762a1bSJed Brown        u(t,0) = 0, u(t,1) = 0,
15c4762a1bSJed Brown    and the initial condition
16c4762a1bSJed Brown        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
17c4762a1bSJed Brown    This is a linear, second-order, parabolic equation.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown    We discretize the right-hand side using finite differences with
20c4762a1bSJed Brown    uniform grid spacing h:
21c4762a1bSJed Brown        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
22c4762a1bSJed Brown    We then demonstrate time evolution using the various TS methods by
23c4762a1bSJed Brown    running the program via
24ab00516aSJared Frazier        mpiexec -n <procs> ./ex4 -ts_type <timestepping solver>
25c4762a1bSJed Brown 
26c4762a1bSJed Brown    We compare the approximate solution with the exact solution, given by
27c4762a1bSJed Brown        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
28c4762a1bSJed Brown                       3*exp(-4*pi*pi*t) * sin(2*pi*x)
29c4762a1bSJed Brown 
30c4762a1bSJed Brown    Notes:
31c4762a1bSJed Brown    This code demonstrates the TS solver interface to two variants of
32c4762a1bSJed Brown    linear problems, u_t = f(u,t), namely
33c4762a1bSJed Brown      - time-dependent f:   f(u,t) is a function of t
34c4762a1bSJed Brown      - time-independent f: f(u,t) is simply f(u)
35c4762a1bSJed Brown 
36c4762a1bSJed Brown     The uniprocessor version of this code is ts/tutorials/ex3.c
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   ------------------------------------------------------------------------- */
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /*
41c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs) to manage
42c4762a1bSJed Brown    the parallel grid.  Include "petscts.h" so that we can use TS solvers.
43c4762a1bSJed Brown    Note that this file 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 <petscdm.h>
52c4762a1bSJed Brown #include <petscdmda.h>
53c4762a1bSJed Brown #include <petscts.h>
54c4762a1bSJed Brown #include <petscdraw.h>
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*
57c4762a1bSJed Brown    User-defined application context - contains data needed by the
58c4762a1bSJed Brown    application-provided call-back routines.
59c4762a1bSJed Brown */
60c4762a1bSJed Brown typedef struct {
61c4762a1bSJed Brown   MPI_Comm    comm;             /* communicator */
62c4762a1bSJed Brown   DM          da;               /* distributed array data structure */
63c4762a1bSJed Brown   Vec         localwork;        /* local ghosted work vector */
64c4762a1bSJed Brown   Vec         u_local;          /* local ghosted approximate solution vector */
65c4762a1bSJed Brown   Vec         solution;         /* global exact solution vector */
66c4762a1bSJed Brown   PetscInt    m;                /* total number of grid points */
67c4762a1bSJed Brown   PetscReal   h;                /* mesh width h = 1/(m-1) */
68c4762a1bSJed Brown   PetscBool   debug;            /* flag (1 indicates activation of debugging printouts) */
69c4762a1bSJed Brown   PetscViewer viewer1, viewer2; /* viewers for the solution and error */
70c4762a1bSJed Brown   PetscReal   norm_2, norm_max; /* error norms */
71c4762a1bSJed Brown } AppCtx;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*
74c4762a1bSJed Brown    User-defined routines
75c4762a1bSJed Brown */
76c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
77c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *);
78c4762a1bSJed Brown extern PetscErrorCode RHSFunctionHeat(TS, PetscReal, Vec, Vec, void *);
79c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
80c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *);
81c4762a1bSJed Brown 
main(int argc,char ** argv)82d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
83d71ae5a4SJacob Faibussowitsch {
84c4762a1bSJed Brown   AppCtx        appctx;               /* user-defined application context */
85c4762a1bSJed Brown   TS            ts;                   /* timestepping context */
86c4762a1bSJed Brown   Mat           A;                    /* matrix data structure */
87c4762a1bSJed Brown   Vec           u;                    /* approximate solution vector */
88c4762a1bSJed Brown   PetscReal     time_total_max = 1.0; /* default max total time */
89c4762a1bSJed Brown   PetscInt      time_steps_max = 100; /* default max timesteps */
90c4762a1bSJed Brown   PetscDraw     draw;                 /* drawing context */
91c4762a1bSJed Brown   PetscInt      steps, m;
92c4762a1bSJed Brown   PetscMPIInt   size;
93c4762a1bSJed Brown   PetscReal     dt, ftime;
94c4762a1bSJed Brown   PetscBool     flg;
95c4762a1bSJed Brown   TSProblemType tsproblem = TS_LINEAR;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown      Initialize program and set problem parameters
99c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100c4762a1bSJed Brown 
101327415f7SBarry Smith   PetscFunctionBeginUser;
102c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
103c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   m = 60;
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
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   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1149566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solving a linear TS problem, number of processors = %d\n", size));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Create vector data structures
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4762a1bSJed Brown   /*
120c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
121c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
122c4762a1bSJed Brown      total grid values spread equally among all the processors.
123c4762a1bSJed Brown   */
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, m, 1, 1, NULL, &appctx.da));
1269566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(appctx.da));
1279566063dSJacob Faibussowitsch   PetscCall(DMSetUp(appctx.da));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /*
130c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
131c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
132c4762a1bSJed Brown      have the same types.
133c4762a1bSJed Brown   */
1349566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(appctx.da, &u));
1359566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(appctx.da, &appctx.u_local));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
139c4762a1bSJed Brown      create global work vector for storing exact solution.
140c4762a1bSJed Brown   */
1419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(appctx.u_local, &appctx.localwork));
1429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
146c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "", 80, 380, 400, 160, &appctx.viewer1));
1499566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer1, 0, &draw));
1509566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
1519566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "", 80, 0, 400, 160, &appctx.viewer2));
1529566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(appctx.viewer2, 0, &draw));
1539566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156c4762a1bSJed Brown      Create timestepping solver context
157c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158c4762a1bSJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   flg = PETSC_FALSE;
1629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-nonlinear", &flg, NULL));
1639566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, flg ? TS_NONLINEAR : TS_LINEAR));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166c4762a1bSJed Brown      Set optional user-defined monitoring routine
167c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1689566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown 
172c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
173c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1769566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m));
1779566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1789566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   flg = PETSC_FALSE;
1819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-time_dependent_rhs", &flg, NULL));
182c4762a1bSJed Brown   if (flg) {
183c4762a1bSJed Brown     /*
184c4762a1bSJed Brown        For linear problems with a time-dependent f(u,t) in the equation
185dd8e379bSPierre Jolivet        u_t = f(u,t), the user provides the discretized right-hand side
186c4762a1bSJed Brown        as a time-dependent matrix.
187c4762a1bSJed Brown     */
1889566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1899566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, A, A, RHSMatrixHeat, &appctx));
190c4762a1bSJed Brown   } else {
191c4762a1bSJed Brown     /*
192c4762a1bSJed Brown        For linear problems with a time-independent f(u) in the equation
193dd8e379bSPierre Jolivet        u_t = f(u), the user provides the discretized right-hand side
194c4762a1bSJed Brown        as a matrix only once, and then sets a null matrix evaluation
195c4762a1bSJed Brown        routine.
196c4762a1bSJed Brown     */
1979566063dSJacob Faibussowitsch     PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx));
1989566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
1999566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, &appctx));
200c4762a1bSJed Brown   }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   if (tsproblem == TS_NONLINEAR) {
203c4762a1bSJed Brown     SNES snes;
2049566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionHeat, &appctx));
2059566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
2069566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESComputeJacobianDefault, NULL));
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Set solution vector and initial timestep
211c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   dt = appctx.h * appctx.h / 2.0;
2149566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
2159566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, u));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218c4762a1bSJed Brown      Customize timestepping solver:
219c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
220c4762a1bSJed Brown        - Set timestepping duration info
221c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
222c4762a1bSJed Brown      For example,
223c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
224c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
225c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226c4762a1bSJed Brown 
2279566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
2289566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
2299566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2309566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233c4762a1bSJed Brown      Solve the problem
234c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /*
237c4762a1bSJed Brown      Evaluate initial conditions
238c4762a1bSJed Brown   */
2399566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /*
242c4762a1bSJed Brown      Run the timestepping solver
243c4762a1bSJed Brown   */
2449566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
2459566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2469566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249c4762a1bSJed Brown      View timestepping solver info
250c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
25163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total timesteps %" PetscInt_FMT ", Final time %g\n", steps, (double)ftime));
2529566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Avg. error (2 norm) = %g Avg. error (max norm) = %g\n", (double)(appctx.norm_2 / steps), (double)(appctx.norm_max / steps)));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
255c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
256c4762a1bSJed Brown      are no longer needed.
257c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2629566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer1));
2639566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&appctx.viewer2));
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.localwork));
2659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
2669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.u_local));
2679566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&appctx.da));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /*
270c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
271c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
272c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
273c4762a1bSJed Brown          options are chosen (e.g., -log_view).
274c4762a1bSJed Brown   */
2759566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
276b122ec5aSJacob Faibussowitsch   return 0;
277c4762a1bSJed Brown }
278c4762a1bSJed Brown /* --------------------------------------------------------------------- */
279c4762a1bSJed Brown /*
280c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
281c4762a1bSJed Brown 
282c4762a1bSJed Brown    Input Parameter:
283c4762a1bSJed Brown    u - uninitialized solution vector (global)
284c4762a1bSJed Brown    appctx - user-defined application context
285c4762a1bSJed Brown 
286c4762a1bSJed Brown    Output Parameter:
287c4762a1bSJed Brown    u - vector with solution at initial time (global)
288c4762a1bSJed Brown */
InitialConditions(Vec u,AppCtx * appctx)289d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
290d71ae5a4SJacob Faibussowitsch {
291c4762a1bSJed Brown   PetscScalar *u_localptr, h = appctx->h;
292c4762a1bSJed Brown   PetscInt     i, mybase, myend;
293c4762a1bSJed Brown 
2943ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
295c4762a1bSJed Brown   /*
296c4762a1bSJed Brown      Determine starting point of each processor's range of
297c4762a1bSJed Brown      grid values.
298c4762a1bSJed Brown   */
2999566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(u, &mybase, &myend));
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   /*
302c4762a1bSJed Brown     Get a pointer to vector data.
303c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
304c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
305c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
306c4762a1bSJed Brown       the array.
307c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
308c4762a1bSJed Brown       C version.  See the users manual for details.
309c4762a1bSJed Brown   */
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &u_localptr));
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /*
313c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
314c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
315c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
316c4762a1bSJed Brown   */
317c4762a1bSJed Brown   for (i = mybase; i < myend; i++) u_localptr[i - mybase] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h);
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   /*
320c4762a1bSJed Brown      Restore vector
321c4762a1bSJed Brown   */
3229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &u_localptr));
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   /*
325c4762a1bSJed Brown      Print debugging information if desired
326c4762a1bSJed Brown   */
327c4762a1bSJed Brown   if (appctx->debug) {
3289566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "initial guess vector\n"));
3299566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
330c4762a1bSJed Brown   }
3313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
332c4762a1bSJed Brown }
333c4762a1bSJed Brown /* --------------------------------------------------------------------- */
334c4762a1bSJed Brown /*
335c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
336c4762a1bSJed Brown 
337c4762a1bSJed Brown    Input Parameters:
338c4762a1bSJed Brown    t - current time
339c4762a1bSJed Brown    solution - vector in which exact solution will be computed
340c4762a1bSJed Brown    appctx - user-defined application context
341c4762a1bSJed Brown 
342c4762a1bSJed Brown    Output Parameter:
343c4762a1bSJed Brown    solution - vector with the newly computed exact solution
344c4762a1bSJed Brown */
ExactSolution(PetscReal t,Vec solution,AppCtx * appctx)345d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx)
346d71ae5a4SJacob Faibussowitsch {
347c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
348c4762a1bSJed Brown   PetscInt     i, mybase, myend;
349c4762a1bSJed Brown 
3503ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
351c4762a1bSJed Brown   /*
352c4762a1bSJed Brown      Determine starting and ending points of each processor's
353c4762a1bSJed Brown      range of grid values
354c4762a1bSJed Brown   */
3559566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(solution, &mybase, &myend));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /*
358c4762a1bSJed Brown      Get a pointer to vector data.
359c4762a1bSJed Brown   */
3609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(solution, &s_localptr));
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   /*
363c4762a1bSJed Brown      Simply write the solution directly into the array locations.
364c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
365c4762a1bSJed Brown   */
3669371c9d4SSatish Balay   ex1 = PetscExpReal(-36. * PETSC_PI * PETSC_PI * t);
3679371c9d4SSatish Balay   ex2 = PetscExpReal(-4. * PETSC_PI * PETSC_PI * t);
3689371c9d4SSatish Balay   sc1 = PETSC_PI * 6. * h;
3699371c9d4SSatish Balay   sc2 = PETSC_PI * 2. * h;
370c4762a1bSJed Brown   for (i = mybase; i < myend; i++) s_localptr[i - mybase] = PetscSinScalar(sc1 * (PetscReal)i) * ex1 + 3. * PetscSinScalar(sc2 * (PetscReal)i) * ex2;
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /*
373c4762a1bSJed Brown      Restore vector
374c4762a1bSJed Brown   */
3759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(solution, &s_localptr));
3763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
377c4762a1bSJed Brown }
378c4762a1bSJed Brown /* --------------------------------------------------------------------- */
379c4762a1bSJed Brown /*
380c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
381c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
382c4762a1bSJed Brown    error in two different norms.
383c4762a1bSJed Brown 
384c4762a1bSJed Brown    Input Parameters:
385c4762a1bSJed Brown    ts     - the timestep context
386c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
387c4762a1bSJed Brown              initial condition)
388c4762a1bSJed Brown    time   - the current time
389c4762a1bSJed Brown    u      - the solution at this timestep
390c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
391c4762a1bSJed Brown             In this case we use the application context which contains
392c4762a1bSJed Brown             information about the problem size, workspace and the exact
393c4762a1bSJed Brown             solution.
394c4762a1bSJed Brown */
Monitor(TS ts,PetscInt step,PetscReal time,Vec u,PetscCtx ctx)395*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, PetscCtx ctx)
396d71ae5a4SJacob Faibussowitsch {
397c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
398c4762a1bSJed Brown   PetscReal norm_2, norm_max;
399c4762a1bSJed Brown 
4003ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
401c4762a1bSJed Brown   /*
402c4762a1bSJed Brown      View a graph of the current iterate
403c4762a1bSJed Brown   */
4049566063dSJacob Faibussowitsch   PetscCall(VecView(u, appctx->viewer2));
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   /*
407c4762a1bSJed Brown      Compute the exact solution
408c4762a1bSJed Brown   */
4099566063dSJacob Faibussowitsch   PetscCall(ExactSolution(time, appctx->solution, appctx));
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   /*
412c4762a1bSJed Brown      Print debugging information if desired
413c4762a1bSJed Brown   */
414c4762a1bSJed Brown   if (appctx->debug) {
4159566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Computed solution vector\n"));
4169566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
4179566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Exact solution vector\n"));
4189566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD));
419c4762a1bSJed Brown   }
420c4762a1bSJed Brown 
421c4762a1bSJed Brown   /*
422c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
423c4762a1bSJed Brown   */
4249566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
4259566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2));
426c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h) * norm_2;
4279566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max));
428c4762a1bSJed Brown   if (norm_2 < 1e-14) norm_2 = 0;
429c4762a1bSJed Brown   if (norm_max < 1e-14) norm_max = 0;
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   /*
432c4762a1bSJed Brown      PetscPrintf() causes only the first processor in this
433c4762a1bSJed Brown      communicator to print the timestep information.
434c4762a1bSJed Brown   */
43563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(appctx->comm, "Timestep %" PetscInt_FMT ": time = %g 2-norm error = %g max norm error = %g\n", step, (double)time, (double)norm_2, (double)norm_max));
436c4762a1bSJed Brown   appctx->norm_2 += norm_2;
437c4762a1bSJed Brown   appctx->norm_max += norm_max;
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /*
440c4762a1bSJed Brown      View a graph of the error
441c4762a1bSJed Brown   */
4429566063dSJacob Faibussowitsch   PetscCall(VecView(appctx->solution, appctx->viewer1));
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   /*
445c4762a1bSJed Brown      Print debugging information if desired
446c4762a1bSJed Brown   */
447c4762a1bSJed Brown   if (appctx->debug) {
4489566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Error vector\n"));
4499566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD));
450c4762a1bSJed Brown   }
4513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
452c4762a1bSJed Brown }
453c4762a1bSJed Brown 
454c4762a1bSJed Brown /* --------------------------------------------------------------------- */
455c4762a1bSJed Brown /*
456c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
457c4762a1bSJed Brown    matrix for the heat equation.
458c4762a1bSJed Brown 
459c4762a1bSJed Brown    Input Parameters:
460c4762a1bSJed Brown    ts - the TS context
461c4762a1bSJed Brown    t - current time
462c4762a1bSJed Brown    global_in - global input vector
463c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
464c4762a1bSJed Brown 
465c4762a1bSJed Brown    Output Parameters:
466c4762a1bSJed Brown    AA - Jacobian matrix
4677addb90fSBarry Smith    BB - optionally different matrix used to construct the preconditioner
468c4762a1bSJed Brown 
469c4762a1bSJed Brown   Notes:
470c4762a1bSJed Brown   RHSMatrixHeat computes entries for the locally owned part of the system.
471c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
472c4762a1bSJed Brown      contiguous chunks of rows across the processors.
473c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
474c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
475c4762a1bSJed Brown      appropriate processor during matrix assembly).
476c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
477c4762a1bSJed Brown      using MatSetValues(); we could alternatively use MatSetValuesLocal().
478c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
479c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
480c4762a1bSJed Brown      in Fortran as well as in C.
481c4762a1bSJed Brown */
RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,PetscCtx ctx)482*2a8381b2SBarry Smith PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec X, Mat AA, Mat BB, PetscCtx ctx)
483d71ae5a4SJacob Faibussowitsch {
484c4762a1bSJed Brown   Mat         A      = AA;            /* Jacobian matrix */
485c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
486c4762a1bSJed Brown   PetscInt    i, mstart, mend, idx[3];
487c4762a1bSJed Brown   PetscScalar v[3], stwo = -2. / (appctx->h * appctx->h), sone = -.5 * stwo;
488c4762a1bSJed Brown 
4893ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
490c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
491c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
492c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
493c4762a1bSJed Brown 
4949566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &mstart, &mend));
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   /*
497c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
498c4762a1bSJed Brown   */
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   if (mstart == 0) { /* first processor only */
501c4762a1bSJed Brown     v[0] = 1.0;
5029566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
503c4762a1bSJed Brown     mstart++;
504c4762a1bSJed Brown   }
505c4762a1bSJed Brown 
506c4762a1bSJed Brown   if (mend == appctx->m) { /* last processor only */
507c4762a1bSJed Brown     mend--;
508c4762a1bSJed Brown     v[0] = 1.0;
5099566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &mend, 1, &mend, v, INSERT_VALUES));
510c4762a1bSJed Brown   }
511c4762a1bSJed Brown 
512c4762a1bSJed Brown   /*
513c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
514c4762a1bSJed Brown      matrix one row at a time.
515c4762a1bSJed Brown   */
5169371c9d4SSatish Balay   v[0] = sone;
5179371c9d4SSatish Balay   v[1] = stwo;
5189371c9d4SSatish Balay   v[2] = sone;
519c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
5209371c9d4SSatish Balay     idx[0] = i - 1;
5219371c9d4SSatish Balay     idx[1] = i;
5229371c9d4SSatish Balay     idx[2] = i + 1;
5239566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
524c4762a1bSJed Brown   }
525c4762a1bSJed Brown 
526c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
527c4762a1bSJed Brown      Complete the matrix assembly process and set some options
528c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
529c4762a1bSJed Brown   /*
530c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
531c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
532c4762a1bSJed Brown      Computations can be done while messages are in transition
533c4762a1bSJed Brown      by placing code between these two statements.
534c4762a1bSJed Brown   */
5359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
537c4762a1bSJed Brown 
538c4762a1bSJed Brown   /*
539c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
540c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
541c4762a1bSJed Brown   */
5429566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
544c4762a1bSJed Brown }
545c4762a1bSJed Brown 
RHSFunctionHeat(TS ts,PetscReal t,Vec globalin,Vec globalout,PetscCtx ctx)546*2a8381b2SBarry Smith PetscErrorCode RHSFunctionHeat(TS ts, PetscReal t, Vec globalin, Vec globalout, PetscCtx ctx)
547d71ae5a4SJacob Faibussowitsch {
548c4762a1bSJed Brown   Mat A;
549c4762a1bSJed Brown 
550c4762a1bSJed Brown   PetscFunctionBeginUser;
5519566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts, &A, NULL, NULL, &ctx));
5529566063dSJacob Faibussowitsch   PetscCall(RHSMatrixHeat(ts, t, globalin, A, NULL, ctx));
5539566063dSJacob Faibussowitsch   /* PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); */
5549566063dSJacob Faibussowitsch   PetscCall(MatMult(A, globalin, globalout));
5553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
556c4762a1bSJed Brown }
557c4762a1bSJed Brown 
558c4762a1bSJed Brown /*TEST
559c4762a1bSJed Brown 
560c4762a1bSJed Brown     test:
561c4762a1bSJed Brown       args: -ts_view -nox
562c4762a1bSJed Brown 
563c4762a1bSJed Brown     test:
564c4762a1bSJed Brown       suffix: 2
565c4762a1bSJed Brown       args: -ts_view -nox
566c4762a1bSJed Brown       nsize: 3
567c4762a1bSJed Brown 
568c4762a1bSJed Brown     test:
569c4762a1bSJed Brown       suffix: 3
570c4762a1bSJed Brown       args: -ts_view -nox -nonlinear
571c4762a1bSJed Brown 
572c4762a1bSJed Brown     test:
573c4762a1bSJed Brown       suffix: 4
574c4762a1bSJed Brown       args: -ts_view -nox -nonlinear
575c4762a1bSJed Brown       nsize: 3
576c4762a1bSJed Brown       timeoutfactor: 3
577c4762a1bSJed Brown 
578c4762a1bSJed Brown     test:
579c4762a1bSJed Brown       suffix: sundials
580e808b789SPatrick Sanan       requires: sundials2
581c4762a1bSJed Brown       args: -nox -ts_type sundials -ts_max_steps 5 -nonlinear
582c4762a1bSJed Brown       nsize: 4
583c4762a1bSJed Brown 
5847324063eSPatrick Sanan     test:
5857324063eSPatrick Sanan       suffix: sundials_dense
5867324063eSPatrick Sanan       requires: sundials2
5877324063eSPatrick Sanan       args: -nox -ts_type sundials -ts_sundials_use_dense -ts_max_steps 5 -nonlinear
5887324063eSPatrick Sanan       nsize: 1
5897324063eSPatrick Sanan 
590c4762a1bSJed Brown TEST*/
591