xref: /petsc/src/ts/tutorials/ex6.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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    Concepts: TS^time-dependent linear problems
11c4762a1bSJed Brown    Concepts: TS^heat equation
12c4762a1bSJed Brown    Concepts: TS^diffusion equation
13c4762a1bSJed Brown    Routines: TSCreate(); TSSetSolution(); TSSetRHSJacobian(), TSSetIJacobian();
14c4762a1bSJed Brown    Routines: TSSetTimeStep(); TSSetMaxTime(); TSMonitorSet();
15c4762a1bSJed Brown    Routines: TSSetFromOptions(); TSStep(); TSDestroy();
16c4762a1bSJed Brown    Routines: TSSetTimeStep(); TSGetTimeStep();
17c4762a1bSJed Brown    Processors: 1
18c4762a1bSJed Brown */
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /* ------------------------------------------------------------------------
21c4762a1bSJed Brown 
22c4762a1bSJed Brown    This program solves the one-dimensional heat equation (also called the
23c4762a1bSJed Brown    diffusion equation),
24c4762a1bSJed Brown        u_t = u_xx,
25c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
26c4762a1bSJed Brown        u(t,0) = 0, u(t,1) = 0,
27c4762a1bSJed Brown    and the initial condition
28c4762a1bSJed Brown        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
29c4762a1bSJed Brown    This is a linear, second-order, parabolic equation.
30c4762a1bSJed Brown 
31c4762a1bSJed Brown    We discretize the right-hand side using finite differences with
32c4762a1bSJed Brown    uniform grid spacing h:
33c4762a1bSJed Brown        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
34c4762a1bSJed Brown    We then demonstrate time evolution using the various TS methods by
35c4762a1bSJed Brown    running the program via
36c4762a1bSJed Brown        ex3 -ts_type <timestepping solver>
37c4762a1bSJed Brown 
38c4762a1bSJed Brown    We compare the approximate solution with the exact solution, given by
39c4762a1bSJed Brown        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
40c4762a1bSJed Brown                       3*exp(-4*pi*pi*t) * sin(2*pi*x)
41c4762a1bSJed Brown 
42c4762a1bSJed Brown    Notes:
43c4762a1bSJed Brown    This code demonstrates the TS solver interface to two variants of
44c4762a1bSJed Brown    linear problems, u_t = f(u,t), namely
45c4762a1bSJed Brown      - time-dependent f:   f(u,t) is a function of t
46c4762a1bSJed Brown      - time-independent f: f(u,t) is simply f(u)
47c4762a1bSJed Brown 
48c4762a1bSJed Brown     The parallel version of this code is ts/tutorials/ex4.c
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   ------------------------------------------------------------------------- */
51c4762a1bSJed Brown 
52c4762a1bSJed Brown /*
53c4762a1bSJed Brown    Include "ts.h" so that we can use TS solvers.  Note that this file
54c4762a1bSJed Brown    automatically includes:
55c4762a1bSJed Brown      petscsys.h  - base PETSc routines   vec.h  - vectors
56c4762a1bSJed Brown      sys.h    - system routines       mat.h  - matrices
57c4762a1bSJed Brown      is.h     - index sets            ksp.h  - Krylov subspace methods
58c4762a1bSJed Brown      viewer.h - viewers               pc.h   - preconditioners
59c4762a1bSJed Brown      snes.h - nonlinear solvers
60c4762a1bSJed Brown */
61c4762a1bSJed Brown 
62c4762a1bSJed Brown #include <petscts.h>
63c4762a1bSJed Brown #include <petscdraw.h>
64c4762a1bSJed Brown 
65c4762a1bSJed Brown /*
66c4762a1bSJed Brown    User-defined application context - contains data needed by the
67c4762a1bSJed Brown    application-provided call-back routines.
68c4762a1bSJed Brown */
69c4762a1bSJed Brown typedef struct {
70c4762a1bSJed Brown   Vec         solution;          /* global exact solution vector */
71c4762a1bSJed Brown   PetscInt    m;                 /* total number of grid points */
72c4762a1bSJed Brown   PetscReal   h;                 /* mesh width h = 1/(m-1) */
73c4762a1bSJed Brown   PetscBool   debug;             /* flag (1 indicates activation of debugging printouts) */
74c4762a1bSJed Brown   PetscViewer viewer1, viewer2;  /* viewers for the solution and error */
75c4762a1bSJed Brown   PetscReal   norm_2, norm_max;  /* error norms */
76c4762a1bSJed Brown } AppCtx;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown /*
79c4762a1bSJed Brown    User-defined routines
80c4762a1bSJed Brown */
81c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
82c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
83c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
84c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
85c4762a1bSJed Brown extern PetscErrorCode MyBCRoutine(TS,PetscReal,Vec,void*);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown int main(int argc,char **argv)
88c4762a1bSJed Brown {
89c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
90c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
91c4762a1bSJed Brown   Mat            A;                      /* matrix data structure */
92c4762a1bSJed Brown   Vec            u;                      /* approximate solution vector */
93c4762a1bSJed Brown   PetscReal      time_total_max = 100.0; /* default max total time */
94c4762a1bSJed Brown   PetscInt       time_steps_max = 100;   /* default max timesteps */
95c4762a1bSJed Brown   PetscDraw      draw;                   /* drawing context */
96c4762a1bSJed Brown   PetscErrorCode ierr;
97c4762a1bSJed Brown   PetscInt       steps, m;
98c4762a1bSJed Brown   PetscMPIInt    size;
99c4762a1bSJed Brown   PetscReal      dt;
100c4762a1bSJed Brown   PetscReal      ftime;
101c4762a1bSJed Brown   PetscBool      flg;
102c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103c4762a1bSJed Brown      Initialize program and set problem parameters
104c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
107c4762a1bSJed Brown   MPI_Comm_size(PETSC_COMM_WORLD,&size);
1083c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   m    = 60;
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   appctx.m        = m;
115c4762a1bSJed Brown   appctx.h        = 1.0/(m-1.0);
116c4762a1bSJed Brown   appctx.norm_2   = 0.0;
117c4762a1bSJed Brown   appctx.norm_max = 0.0;
118c4762a1bSJed Brown 
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n"));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4762a1bSJed Brown      Create vector data structures
123c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /*
126c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
127c4762a1bSJed Brown   */
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&u));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&appctx.solution));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
133c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134c4762a1bSJed Brown 
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1));
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer1,0,&draw));
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetDoubleBuffer(draw));
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer2,0,&draw));
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetDoubleBuffer(draw));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143c4762a1bSJed Brown      Create timestepping solver context
144c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145c4762a1bSJed Brown 
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_LINEAR));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown      Set optional user-defined monitoring routine
151c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152c4762a1bSJed Brown 
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSMonitorSet(ts,Monitor,&appctx,NULL));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156c4762a1bSJed Brown 
157c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
158c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159c4762a1bSJed Brown 
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
164c4762a1bSJed Brown 
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-time_dependent_rhs",&flg));
166c4762a1bSJed Brown   if (flg) {
167c4762a1bSJed Brown     /*
168c4762a1bSJed Brown        For linear problems with a time-dependent f(u,t) in the equation
169c4762a1bSJed Brown        u_t = f(u,t), the user provides the discretized right-hand-side
170c4762a1bSJed Brown        as a time-dependent matrix.
171c4762a1bSJed Brown     */
172*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx));
173*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx));
174c4762a1bSJed Brown   } else {
175c4762a1bSJed Brown     /*
176c4762a1bSJed Brown        For linear problems with a time-independent f(u) in the equation
177c4762a1bSJed Brown        u_t = f(u), the user provides the discretized right-hand-side
178c4762a1bSJed Brown        as a matrix only once, and then sets a null matrix evaluation
179c4762a1bSJed Brown        routine.
180c4762a1bSJed Brown     */
181*5f80ce2aSJacob Faibussowitsch     CHKERRQ(RHSMatrixHeat(ts,0.0,u,A,A,&appctx));
182*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx));
183*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx));
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown      Set solution vector and initial timestep
188c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   dt   = appctx.h*appctx.h/2.0;
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,u));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195c4762a1bSJed Brown      Customize timestepping solver:
196c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
197c4762a1bSJed Brown        - Set timestepping duration info
198c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
199c4762a1bSJed Brown      For example,
200c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
201c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
202c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203c4762a1bSJed Brown 
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,time_steps_max));
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,time_total_max));
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Solve the problem
211c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /*
214c4762a1bSJed Brown      Evaluate initial conditions
215c4762a1bSJed Brown   */
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InitialConditions(u,&appctx));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /*
219c4762a1bSJed Brown      Run the timestepping solver
220c4762a1bSJed Brown   */
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,u));
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226c4762a1bSJed Brown      View timestepping solver info
227c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228c4762a1bSJed Brown 
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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)));
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSView(ts,PETSC_VIEWER_STDOUT_SELF));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
234c4762a1bSJed Brown      are no longer needed.
235c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236c4762a1bSJed Brown 
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
238*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&appctx.viewer1));
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&appctx.viewer2));
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&appctx.solution));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /*
245c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
246c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
247c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
248c4762a1bSJed Brown          options are chosen (e.g., -log_view).
249c4762a1bSJed Brown   */
250c4762a1bSJed Brown   ierr = PetscFinalize();
251c4762a1bSJed Brown   return ierr;
252c4762a1bSJed Brown }
253c4762a1bSJed Brown /* --------------------------------------------------------------------- */
254c4762a1bSJed Brown /*
255c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
256c4762a1bSJed Brown 
257c4762a1bSJed Brown    Input Parameter:
258c4762a1bSJed Brown    u - uninitialized solution vector (global)
259c4762a1bSJed Brown    appctx - user-defined application context
260c4762a1bSJed Brown 
261c4762a1bSJed Brown    Output Parameter:
262c4762a1bSJed Brown    u - vector with solution at initial time (global)
263c4762a1bSJed Brown */
264c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
265c4762a1bSJed Brown {
266c4762a1bSJed Brown   PetscScalar    *u_localptr;
267c4762a1bSJed Brown   PetscInt       i;
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /*
270c4762a1bSJed Brown     Get a pointer to vector data.
271c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
272c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
273c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
274c4762a1bSJed Brown       the array.
275c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
276c4762a1bSJed Brown       C version.  See the users manual for details.
277c4762a1bSJed Brown   */
278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(u,&u_localptr));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   /*
281c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
282c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
283c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
284c4762a1bSJed Brown   */
285c4762a1bSJed 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);
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /*
288c4762a1bSJed Brown      Restore vector
289c4762a1bSJed Brown   */
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(u,&u_localptr));
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   /*
293c4762a1bSJed Brown      Print debugging information if desired
294c4762a1bSJed Brown   */
295c4762a1bSJed Brown   if (appctx->debug) {
296*5f80ce2aSJacob Faibussowitsch      CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF));
297c4762a1bSJed Brown   }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   return 0;
300c4762a1bSJed Brown }
301c4762a1bSJed Brown /* --------------------------------------------------------------------- */
302c4762a1bSJed Brown /*
303c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
304c4762a1bSJed Brown 
305c4762a1bSJed Brown    Input Parameters:
306c4762a1bSJed Brown    t - current time
307c4762a1bSJed Brown    solution - vector in which exact solution will be computed
308c4762a1bSJed Brown    appctx - user-defined application context
309c4762a1bSJed Brown 
310c4762a1bSJed Brown    Output Parameter:
311c4762a1bSJed Brown    solution - vector with the newly computed exact solution
312c4762a1bSJed Brown */
313c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
314c4762a1bSJed Brown {
315c4762a1bSJed Brown   PetscScalar    *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
316c4762a1bSJed Brown   PetscInt       i;
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /*
319c4762a1bSJed Brown      Get a pointer to vector data.
320c4762a1bSJed Brown   */
321*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(solution,&s_localptr));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /*
324c4762a1bSJed Brown      Simply write the solution directly into the array locations.
325c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
326c4762a1bSJed Brown   */
327c4762a1bSJed Brown   ex1 = PetscExpReal(-36.*PETSC_PI*PETSC_PI*t); ex2 = PetscExpReal(-4.*PETSC_PI*PETSC_PI*t);
328c4762a1bSJed Brown   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
329c4762a1bSJed Brown   for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinReal(PetscRealPart(sc1)*(PetscReal)i)*ex1 + 3.*PetscSinReal(PetscRealPart(sc2)*(PetscReal)i)*ex2;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   /*
332c4762a1bSJed Brown      Restore vector
333c4762a1bSJed Brown   */
334*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(solution,&s_localptr));
335c4762a1bSJed Brown   return 0;
336c4762a1bSJed Brown }
337c4762a1bSJed Brown /* --------------------------------------------------------------------- */
338c4762a1bSJed Brown /*
339c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
340c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
341c4762a1bSJed Brown    error in two different norms.
342c4762a1bSJed Brown 
343c4762a1bSJed Brown    This example also demonstrates changing the timestep via TSSetTimeStep().
344c4762a1bSJed Brown 
345c4762a1bSJed Brown    Input Parameters:
346c4762a1bSJed Brown    ts     - the timestep context
347c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
348c4762a1bSJed Brown              initial condition)
349c4762a1bSJed Brown    crtime  - the current time
350c4762a1bSJed Brown    u      - the solution at this timestep
351c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
352c4762a1bSJed Brown             In this case we use the application context which contains
353c4762a1bSJed Brown             information about the problem size, workspace and the exact
354c4762a1bSJed Brown             solution.
355c4762a1bSJed Brown */
356c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal crtime,Vec u,void *ctx)
357c4762a1bSJed Brown {
358c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
359c4762a1bSJed Brown   PetscReal      norm_2, norm_max, dt, dttol;
360c4762a1bSJed Brown   PetscBool      flg;
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   /*
363c4762a1bSJed Brown      View a graph of the current iterate
364c4762a1bSJed Brown   */
365*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(u,appctx->viewer2));
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   /*
368c4762a1bSJed Brown      Compute the exact solution
369c4762a1bSJed Brown   */
370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ExactSolution(crtime,appctx->solution,appctx));
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /*
373c4762a1bSJed Brown      Print debugging information if desired
374c4762a1bSJed Brown   */
375c4762a1bSJed Brown   if (appctx->debug) {
376*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n"));
377*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF));
378*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n"));
379*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF));
380c4762a1bSJed Brown   }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   /*
383c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
384c4762a1bSJed Brown   */
385*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(appctx->solution,-1.0,u));
386*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(appctx->solution,NORM_2,&norm_2));
387c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
388*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(appctx->solution,NORM_MAX,&norm_max));
389c4762a1bSJed Brown 
390*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
391c4762a1bSJed Brown   if (norm_2 > 1.e-2) {
392*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Timestep %D: 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));
393c4762a1bSJed Brown   }
394c4762a1bSJed Brown   appctx->norm_2   += norm_2;
395c4762a1bSJed Brown   appctx->norm_max += norm_max;
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   dttol = .0001;
398*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,&flg));
399c4762a1bSJed Brown   if (dt < dttol) {
400c4762a1bSJed Brown     dt  *= .999;
401*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetTimeStep(ts,dt));
402c4762a1bSJed Brown   }
403c4762a1bSJed Brown 
404c4762a1bSJed Brown   /*
405c4762a1bSJed Brown      View a graph of the error
406c4762a1bSJed Brown   */
407*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(appctx->solution,appctx->viewer1));
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   /*
410c4762a1bSJed Brown      Print debugging information if desired
411c4762a1bSJed Brown   */
412c4762a1bSJed Brown   if (appctx->debug) {
413*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error vector\n"));
414*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF));
415c4762a1bSJed Brown   }
416c4762a1bSJed Brown 
417c4762a1bSJed Brown   return 0;
418c4762a1bSJed Brown }
419c4762a1bSJed Brown /* --------------------------------------------------------------------- */
420c4762a1bSJed Brown /*
421c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
422c4762a1bSJed Brown    matrix for the heat equation.
423c4762a1bSJed Brown 
424c4762a1bSJed Brown    Input Parameters:
425c4762a1bSJed Brown    ts - the TS context
426c4762a1bSJed Brown    t - current time
427c4762a1bSJed Brown    global_in - global input vector
428c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
429c4762a1bSJed Brown 
430c4762a1bSJed Brown    Output Parameters:
431c4762a1bSJed Brown    AA - Jacobian matrix
432c4762a1bSJed Brown    BB - optionally different preconditioning matrix
433c4762a1bSJed Brown    str - flag indicating matrix structure
434c4762a1bSJed Brown 
435c4762a1bSJed Brown    Notes:
436c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
437c4762a1bSJed Brown    in Fortran as well as in C.
438c4762a1bSJed Brown */
439c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
440c4762a1bSJed Brown {
441c4762a1bSJed Brown   Mat            A       = AA;                /* Jacobian matrix */
442c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*) ctx;      /* user-defined application context */
443c4762a1bSJed Brown   PetscInt       mstart  = 0;
444c4762a1bSJed Brown   PetscInt       mend    = appctx->m;
445c4762a1bSJed Brown   PetscInt       i, idx[3];
446c4762a1bSJed Brown   PetscScalar    v[3], stwo = -2./(appctx->h*appctx->h), sone = -.5*stwo;
447c4762a1bSJed Brown 
448c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
449c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
450c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
451c4762a1bSJed Brown   /*
452c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
453c4762a1bSJed Brown   */
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   mstart = 0;
456c4762a1bSJed Brown   v[0]   = 1.0;
457*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES));
458c4762a1bSJed Brown   mstart++;
459c4762a1bSJed Brown 
460c4762a1bSJed Brown   mend--;
461c4762a1bSJed Brown   v[0] = 1.0;
462*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES));
463c4762a1bSJed Brown 
464c4762a1bSJed Brown   /*
465c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
466c4762a1bSJed Brown      matrix one row at a time.
467c4762a1bSJed Brown   */
468c4762a1bSJed Brown   v[0] = sone; v[1] = stwo; v[2] = sone;
469c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
470c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
471*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES));
472c4762a1bSJed Brown   }
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
475c4762a1bSJed Brown      Complete the matrix assembly process and set some options
476c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
477c4762a1bSJed Brown   /*
478c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
479c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
480c4762a1bSJed Brown      Computations can be done while messages are in transition
481c4762a1bSJed Brown      by placing code between these two statements.
482c4762a1bSJed Brown   */
483*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
484*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
485c4762a1bSJed Brown 
486c4762a1bSJed Brown   /*
487c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
488c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
489c4762a1bSJed Brown   */
490*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   return 0;
493c4762a1bSJed Brown }
494c4762a1bSJed Brown /* --------------------------------------------------------------------- */
495c4762a1bSJed Brown /*
496c4762a1bSJed Brown    Input Parameters:
497c4762a1bSJed Brown    ts - the TS context
498c4762a1bSJed Brown    t - current time
499c4762a1bSJed Brown    f - function
500c4762a1bSJed Brown    ctx - optional user-defined context, as set by TSetBCFunction()
501c4762a1bSJed Brown  */
502c4762a1bSJed Brown PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
503c4762a1bSJed Brown {
504c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*) ctx;      /* user-defined application context */
505c4762a1bSJed Brown   PetscInt       m = appctx->m;
506c4762a1bSJed Brown   PetscScalar    *fa;
507c4762a1bSJed Brown 
508*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(f,&fa));
509c4762a1bSJed Brown   fa[0]   = 0.0;
510c4762a1bSJed Brown   fa[m-1] = 1.0;
511*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(f,&fa));
512*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"t=%g\n",(double)t));
513c4762a1bSJed Brown 
514c4762a1bSJed Brown   return 0;
515c4762a1bSJed Brown }
516c4762a1bSJed Brown 
517c4762a1bSJed Brown /*TEST
518c4762a1bSJed Brown 
519c4762a1bSJed Brown     test:
520c4762a1bSJed Brown       args: -nox -ts_max_steps 4
521c4762a1bSJed Brown 
522c4762a1bSJed Brown TEST*/
523