xref: /petsc/src/ts/tutorials/ex3.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   -use_ifunc          : Use IFunction/IJacobian interface\n\
7c4762a1bSJed Brown   -debug              : Activate debugging printouts\n\
8c4762a1bSJed Brown   -nox                : Deactivate x-window graphics\n\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown    Concepts: TS^time-dependent linear problems
12c4762a1bSJed Brown    Concepts: TS^heat equation
13c4762a1bSJed Brown    Concepts: TS^diffusion equation
14c4762a1bSJed Brown    Processors: 1
15c4762a1bSJed Brown */
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /* ------------------------------------------------------------------------
18c4762a1bSJed Brown 
19c4762a1bSJed Brown    This program solves the one-dimensional heat equation (also called the
20c4762a1bSJed Brown    diffusion equation),
21c4762a1bSJed Brown        u_t = u_xx,
22c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
23c4762a1bSJed Brown        u(t,0) = 0, u(t,1) = 0,
24c4762a1bSJed Brown    and the initial condition
25c4762a1bSJed Brown        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
26c4762a1bSJed Brown    This is a linear, second-order, parabolic equation.
27c4762a1bSJed Brown 
28c4762a1bSJed Brown    We discretize the right-hand side using finite differences with
29c4762a1bSJed Brown    uniform grid spacing h:
30c4762a1bSJed Brown        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
31c4762a1bSJed Brown    We then demonstrate time evolution using the various TS methods by
32c4762a1bSJed Brown    running the program via
33c4762a1bSJed Brown        ex3 -ts_type <timestepping solver>
34c4762a1bSJed Brown 
35c4762a1bSJed Brown    We compare the approximate solution with the exact solution, given by
36c4762a1bSJed Brown        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
37c4762a1bSJed Brown                       3*exp(-4*pi*pi*t) * sin(2*pi*x)
38c4762a1bSJed Brown 
39c4762a1bSJed Brown    Notes:
40c4762a1bSJed Brown    This code demonstrates the TS solver interface to two variants of
41c4762a1bSJed Brown    linear problems, u_t = f(u,t), namely
42c4762a1bSJed Brown      - time-dependent f:   f(u,t) is a function of t
43c4762a1bSJed Brown      - time-independent f: f(u,t) is simply f(u)
44c4762a1bSJed Brown 
45c4762a1bSJed Brown     The parallel version of this code is ts/tutorials/ex4.c
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   ------------------------------------------------------------------------- */
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*
50c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this file
51c4762a1bSJed Brown    automatically includes:
52c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h  - vectors
53c4762a1bSJed Brown      petscmat.h  - matrices
54c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
55c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h   - preconditioners
56c4762a1bSJed Brown      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
57c4762a1bSJed Brown */
58c4762a1bSJed Brown 
59c4762a1bSJed Brown #include <petscts.h>
60c4762a1bSJed Brown #include <petscdraw.h>
61c4762a1bSJed Brown 
62c4762a1bSJed Brown /*
63c4762a1bSJed Brown    User-defined application context - contains data needed by the
64c4762a1bSJed Brown    application-provided call-back routines.
65c4762a1bSJed Brown */
66c4762a1bSJed Brown typedef struct {
67c4762a1bSJed Brown   Vec         solution;          /* global exact solution vector */
68c4762a1bSJed Brown   PetscInt    m;                 /* total number of grid points */
69c4762a1bSJed Brown   PetscReal   h;                 /* mesh width h = 1/(m-1) */
70c4762a1bSJed Brown   PetscBool   debug;             /* flag (1 indicates activation of debugging printouts) */
71c4762a1bSJed Brown   PetscViewer viewer1,viewer2;   /* viewers for the solution and error */
72c4762a1bSJed Brown   PetscReal   norm_2,norm_max;   /* error norms */
73c4762a1bSJed Brown   Mat         A;                 /* RHS mat, used with IFunction interface */
74c4762a1bSJed Brown   PetscReal   oshift;            /* old shift applied, prevent to recompute the IJacobian */
75c4762a1bSJed Brown } AppCtx;
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*
78c4762a1bSJed Brown    User-defined routines
79c4762a1bSJed Brown */
80c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
81c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
82c4762a1bSJed Brown extern PetscErrorCode IFunctionHeat(TS,PetscReal,Vec,Vec,Vec,void*);
83c4762a1bSJed Brown extern PetscErrorCode IJacobianHeat(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
84c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
85c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
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   PetscBool      flg,flg_string;
101c4762a1bSJed Brown 
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;
107*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(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   flg_string = PETSC_FALSE;
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_string_viewer",&flg_string,NULL));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   appctx.m        = m;
117c4762a1bSJed Brown   appctx.h        = 1.0/(m-1.0);
118c4762a1bSJed Brown   appctx.norm_2   = 0.0;
119c4762a1bSJed Brown   appctx.norm_max = 0.0;
120c4762a1bSJed Brown 
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n"));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown      Create vector data structures
125c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /*
128c4762a1bSJed Brown      Create vector data structures for approximate and exact solutions
129c4762a1bSJed Brown   */
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&u));
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&appctx.solution));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
135c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136c4762a1bSJed Brown 
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1));
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer1,0,&draw));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetDoubleBuffer(draw));
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2));
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer2,0,&draw));
142*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetDoubleBuffer(draw));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Create timestepping solver context
146c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147c4762a1bSJed Brown 
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_LINEAR));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152c4762a1bSJed Brown      Set optional user-defined monitoring routine
153c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   if (!flg_string) {
156*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,Monitor,&appctx,NULL));
157c4762a1bSJed Brown   }
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160c4762a1bSJed Brown 
161c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
162c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163c4762a1bSJed Brown 
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m));
166*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
167*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   flg  = PETSC_FALSE;
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_ifunc",&flg,NULL));
171c4762a1bSJed Brown   if (!flg) {
172c4762a1bSJed Brown     appctx.A = NULL;
173*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL));
174c4762a1bSJed Brown     if (flg) {
175c4762a1bSJed Brown       /*
176c4762a1bSJed Brown          For linear problems with a time-dependent f(u,t) in the equation
177c4762a1bSJed Brown          u_t = f(u,t), the user provides the discretized right-hand-side
178c4762a1bSJed Brown          as a time-dependent matrix.
179c4762a1bSJed Brown       */
180*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx));
181*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx));
182c4762a1bSJed Brown     } else {
183c4762a1bSJed Brown       /*
184c4762a1bSJed Brown          For linear problems with a time-independent f(u) in the equation
185c4762a1bSJed Brown          u_t = f(u), the user provides the discretized right-hand-side
186c4762a1bSJed Brown          as a matrix only once, and then sets the special Jacobian evaluation
187c4762a1bSJed Brown          routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian.
188c4762a1bSJed Brown       */
189*5f80ce2aSJacob Faibussowitsch       CHKERRQ(RHSMatrixHeat(ts,0.0,u,A,A,&appctx));
190*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx));
191*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx));
192c4762a1bSJed Brown     }
193c4762a1bSJed Brown   } else {
194c4762a1bSJed Brown     Mat J;
195c4762a1bSJed Brown 
196*5f80ce2aSJacob Faibussowitsch     CHKERRQ(RHSMatrixHeat(ts,0.0,u,A,A,&appctx));
197*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&J));
198*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIFunction(ts,NULL,IFunctionHeat,&appctx));
199*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIJacobian(ts,J,J,IJacobianHeat,&appctx));
200*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&J));
201c4762a1bSJed Brown 
202*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectReference((PetscObject)A));
203c4762a1bSJed Brown     appctx.A = A;
204c4762a1bSJed Brown     appctx.oshift = PETSC_MIN_REAL;
205c4762a1bSJed Brown   }
206c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207c4762a1bSJed Brown      Set solution vector and initial timestep
208c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   dt   = appctx.h*appctx.h/2.0;
211*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214c4762a1bSJed Brown      Customize timestepping solver:
215c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
216c4762a1bSJed Brown        - Set timestepping duration info
217c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
218c4762a1bSJed Brown      For example,
219c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
220c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
221c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
222c4762a1bSJed Brown 
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,time_steps_max));
224*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,time_total_max));
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229c4762a1bSJed Brown      Solve the problem
230c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /*
233c4762a1bSJed Brown      Evaluate initial conditions
234c4762a1bSJed Brown   */
235*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InitialConditions(u,&appctx));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /*
238c4762a1bSJed Brown      Run the timestepping solver
239c4762a1bSJed Brown   */
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,u));
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown      View timestepping solver info
245c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246c4762a1bSJed Brown 
247*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)));
248c4762a1bSJed Brown   if (!flg_string) {
249*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSView(ts,PETSC_VIEWER_STDOUT_SELF));
250c4762a1bSJed Brown   } else {
251c4762a1bSJed Brown     PetscViewer stringviewer;
252c4762a1bSJed Brown     char        string[512];
253c4762a1bSJed Brown     const char  *outstring;
254c4762a1bSJed Brown 
255*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerStringOpen(PETSC_COMM_WORLD,string,sizeof(string),&stringviewer));
256*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSView(ts,stringviewer));
257*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerStringGetStringRead(stringviewer,&outstring,NULL));
2583c633725SBarry Smith     PetscCheck((char*)outstring == (char*)string,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"String returned from viewer does not equal original string");
259*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Output from string viewer:%s\n",outstring));
260*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&stringviewer));
261c4762a1bSJed Brown   }
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
265c4762a1bSJed Brown      are no longer needed.
266c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
267c4762a1bSJed Brown 
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
270*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
271*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&appctx.viewer1));
272*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&appctx.viewer2));
273*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&appctx.solution));
274*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&appctx.A));
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   /*
277c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
278c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
279c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
280c4762a1bSJed Brown          options are chosen (e.g., -log_view).
281c4762a1bSJed Brown   */
282c4762a1bSJed Brown   ierr = PetscFinalize();
283c4762a1bSJed Brown   return ierr;
284c4762a1bSJed Brown }
285c4762a1bSJed Brown /* --------------------------------------------------------------------- */
286c4762a1bSJed Brown /*
287c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
288c4762a1bSJed Brown 
289c4762a1bSJed Brown    Input Parameter:
290c4762a1bSJed Brown    u - uninitialized solution vector (global)
291c4762a1bSJed Brown    appctx - user-defined application context
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    Output Parameter:
294c4762a1bSJed Brown    u - vector with solution at initial time (global)
295c4762a1bSJed Brown */
296c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
297c4762a1bSJed Brown {
298c4762a1bSJed Brown   PetscScalar    *u_localptr,h = appctx->h;
299c4762a1bSJed Brown   PetscInt       i;
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   */
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(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=0; i<appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   /*
320c4762a1bSJed Brown      Restore vector
321c4762a1bSJed Brown   */
322*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(u,&u_localptr));
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   /*
325c4762a1bSJed Brown      Print debugging information if desired
326c4762a1bSJed Brown   */
327c4762a1bSJed Brown   if (appctx->debug) {
328*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Initial guess vector\n"));
329*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF));
330c4762a1bSJed Brown   }
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   return 0;
333c4762a1bSJed Brown }
334c4762a1bSJed Brown /* --------------------------------------------------------------------- */
335c4762a1bSJed Brown /*
336c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
337c4762a1bSJed Brown 
338c4762a1bSJed Brown    Input Parameters:
339c4762a1bSJed Brown    t - current time
340c4762a1bSJed Brown    solution - vector in which exact solution will be computed
341c4762a1bSJed Brown    appctx - user-defined application context
342c4762a1bSJed Brown 
343c4762a1bSJed Brown    Output Parameter:
344c4762a1bSJed Brown    solution - vector with the newly computed exact solution
345c4762a1bSJed Brown */
346c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
347c4762a1bSJed Brown {
348c4762a1bSJed Brown   PetscScalar    *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
349c4762a1bSJed Brown   PetscInt       i;
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   /*
352c4762a1bSJed Brown      Get a pointer to vector data.
353c4762a1bSJed Brown   */
354*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(solution,&s_localptr));
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   /*
357c4762a1bSJed Brown      Simply write the solution directly into the array locations.
358c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
359c4762a1bSJed Brown   */
360c4762a1bSJed Brown   ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc);
361c4762a1bSJed Brown   ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
362c4762a1bSJed Brown   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
363c4762a1bSJed Brown   for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   /*
366c4762a1bSJed Brown      Restore vector
367c4762a1bSJed Brown   */
368*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(solution,&s_localptr));
369c4762a1bSJed Brown   return 0;
370c4762a1bSJed Brown }
371c4762a1bSJed Brown /* --------------------------------------------------------------------- */
372c4762a1bSJed Brown /*
373c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
374c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
375c4762a1bSJed Brown    error in two different norms.
376c4762a1bSJed Brown 
377c4762a1bSJed Brown    This example also demonstrates changing the timestep via TSSetTimeStep().
378c4762a1bSJed Brown 
379c4762a1bSJed Brown    Input Parameters:
380c4762a1bSJed Brown    ts     - the timestep context
381c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
382c4762a1bSJed Brown              initial condition)
383c4762a1bSJed Brown    time   - the current time
384c4762a1bSJed Brown    u      - the solution at this timestep
385c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
386c4762a1bSJed Brown             In this case we use the application context which contains
387c4762a1bSJed Brown             information about the problem size, workspace and the exact
388c4762a1bSJed Brown             solution.
389c4762a1bSJed Brown */
390c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
391c4762a1bSJed Brown {
392c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
393c4762a1bSJed Brown   PetscReal      norm_2,norm_max,dt,dttol;
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   /*
396c4762a1bSJed Brown      View a graph of the current iterate
397c4762a1bSJed Brown   */
398*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(u,appctx->viewer2));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   /*
401c4762a1bSJed Brown      Compute the exact solution
402c4762a1bSJed Brown   */
403*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ExactSolution(time,appctx->solution,appctx));
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   /*
406c4762a1bSJed Brown      Print debugging information if desired
407c4762a1bSJed Brown   */
408c4762a1bSJed Brown   if (appctx->debug) {
409*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n"));
410*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF));
411*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n"));
412*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF));
413c4762a1bSJed Brown   }
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   /*
416c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
417c4762a1bSJed Brown   */
418*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(appctx->solution,-1.0,u));
419*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(appctx->solution,NORM_2,&norm_2));
420c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
421*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(appctx->solution,NORM_MAX,&norm_max));
422c4762a1bSJed Brown 
423*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
424*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Timestep %3D: step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)dt,(double)time,(double)norm_2,(double)norm_max));
425c4762a1bSJed Brown 
426c4762a1bSJed Brown   appctx->norm_2   += norm_2;
427c4762a1bSJed Brown   appctx->norm_max += norm_max;
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   dttol = .0001;
430*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,NULL));
431c4762a1bSJed Brown   if (dt < dttol) {
432c4762a1bSJed Brown     dt  *= .999;
433*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetTimeStep(ts,dt));
434c4762a1bSJed Brown   }
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /*
437c4762a1bSJed Brown      View a graph of the error
438c4762a1bSJed Brown   */
439*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(appctx->solution,appctx->viewer1));
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /*
442c4762a1bSJed Brown      Print debugging information if desired
443c4762a1bSJed Brown   */
444c4762a1bSJed Brown   if (appctx->debug) {
445*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error vector\n"));
446*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF));
447c4762a1bSJed Brown   }
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   return 0;
450c4762a1bSJed Brown }
451c4762a1bSJed Brown /* --------------------------------------------------------------------- */
452c4762a1bSJed Brown /*
453c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
454c4762a1bSJed Brown    matrix for the heat equation.
455c4762a1bSJed Brown 
456c4762a1bSJed Brown    Input Parameters:
457c4762a1bSJed Brown    ts - the TS context
458c4762a1bSJed Brown    t - current time
459c4762a1bSJed Brown    global_in - global input vector
460c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
461c4762a1bSJed Brown 
462c4762a1bSJed Brown    Output Parameters:
463c4762a1bSJed Brown    AA - Jacobian matrix
464c4762a1bSJed Brown    BB - optionally different preconditioning matrix
465c4762a1bSJed Brown    str - flag indicating matrix structure
466c4762a1bSJed Brown 
467c4762a1bSJed Brown    Notes:
468c4762a1bSJed Brown    Recall that MatSetValues() uses 0-based row and column numbers
469c4762a1bSJed Brown    in Fortran as well as in C.
470c4762a1bSJed Brown */
471c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
472c4762a1bSJed Brown {
473c4762a1bSJed Brown   Mat            A       = AA;                /* Jacobian matrix */
474c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
475c4762a1bSJed Brown   PetscInt       mstart  = 0;
476c4762a1bSJed Brown   PetscInt       mend    = appctx->m;
477c4762a1bSJed Brown   PetscInt       i,idx[3];
478c4762a1bSJed Brown   PetscScalar    v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;
479c4762a1bSJed Brown 
480c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
481c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
482c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
483c4762a1bSJed Brown   /*
484c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
485c4762a1bSJed Brown   */
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   mstart = 0;
488c4762a1bSJed Brown   v[0]   = 1.0;
489*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES));
490c4762a1bSJed Brown   mstart++;
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   mend--;
493c4762a1bSJed Brown   v[0] = 1.0;
494*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES));
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   /*
497c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
498c4762a1bSJed Brown      matrix one row at a time.
499c4762a1bSJed Brown   */
500c4762a1bSJed Brown   v[0] = sone; v[1] = stwo; v[2] = sone;
501c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
502c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
503*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES));
504c4762a1bSJed Brown   }
505c4762a1bSJed Brown 
506c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
507c4762a1bSJed Brown      Complete the matrix assembly process and set some options
508c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
509c4762a1bSJed Brown   /*
510c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
511c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
512c4762a1bSJed Brown      Computations can be done while messages are in transition
513c4762a1bSJed Brown      by placing code between these two statements.
514c4762a1bSJed Brown   */
515*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
516*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
517c4762a1bSJed Brown 
518c4762a1bSJed Brown   /*
519c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
520c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
521c4762a1bSJed Brown   */
522*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
523c4762a1bSJed Brown 
524c4762a1bSJed Brown   return 0;
525c4762a1bSJed Brown }
526c4762a1bSJed Brown 
527c4762a1bSJed Brown PetscErrorCode IFunctionHeat(TS ts,PetscReal t,Vec X,Vec Xdot,Vec r,void *ctx)
528c4762a1bSJed Brown {
529c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
530c4762a1bSJed Brown 
531*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(appctx->A,X,r));
532*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAYPX(r,-1.0,Xdot));
533c4762a1bSJed Brown   return 0;
534c4762a1bSJed Brown }
535c4762a1bSJed Brown 
536c4762a1bSJed Brown PetscErrorCode IJacobianHeat(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal s,Mat A,Mat B,void *ctx)
537c4762a1bSJed Brown {
538c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
539c4762a1bSJed Brown 
540c4762a1bSJed Brown   if (appctx->oshift == s) return 0;
541*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCopy(appctx->A,A,SAME_NONZERO_PATTERN));
542*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,-1));
543*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,s));
544*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCopy(A,B,SAME_NONZERO_PATTERN));
545c4762a1bSJed Brown   appctx->oshift = s;
546c4762a1bSJed Brown   return 0;
547c4762a1bSJed Brown }
548c4762a1bSJed Brown 
549c4762a1bSJed Brown /*TEST
550c4762a1bSJed Brown 
551c4762a1bSJed Brown     test:
552c4762a1bSJed Brown       args: -nox -ts_type ssp -ts_dt 0.0005
553c4762a1bSJed Brown 
554c4762a1bSJed Brown     test:
555c4762a1bSJed Brown       suffix: 2
556c4762a1bSJed Brown       args: -nox -ts_type ssp -ts_dt 0.0005 -time_dependent_rhs 1
557c4762a1bSJed Brown 
558c4762a1bSJed Brown     test:
559c4762a1bSJed Brown       suffix: 3
560c4762a1bSJed Brown       args:  -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason
561c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
562c4762a1bSJed Brown       requires: !single
563c4762a1bSJed Brown 
564c4762a1bSJed Brown     test:
565c4762a1bSJed Brown       suffix: 4
566c4762a1bSJed Brown       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason
567c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
568c4762a1bSJed Brown 
569c4762a1bSJed Brown     test:
570c4762a1bSJed Brown       suffix: 5
571c4762a1bSJed Brown       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs
572c4762a1bSJed Brown       filter: sed "s/ATOL/RTOL/g"
573c4762a1bSJed Brown 
574c4762a1bSJed Brown     test:
575c4762a1bSJed Brown       requires: !single
576c4762a1bSJed Brown       suffix: pod_guess
577c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason
578c4762a1bSJed Brown 
579c4762a1bSJed Brown     test:
580c4762a1bSJed Brown       requires: !single
581c4762a1bSJed Brown       suffix: pod_guess_Ainner
582c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -ksp_guess_pod_Ainner -pc_type none -ksp_converged_reason
583c4762a1bSJed Brown 
584c4762a1bSJed Brown     test:
585c4762a1bSJed Brown       requires: !single
586c4762a1bSJed Brown       suffix: fischer_guess
587c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason
588c4762a1bSJed Brown 
589c4762a1bSJed Brown     test:
590c4762a1bSJed Brown       requires: !single
591c4762a1bSJed Brown       suffix: fischer_guess_2
592c4762a1bSJed Brown       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 2,10 -pc_type none -ksp_converged_reason
593c4762a1bSJed Brown 
594c4762a1bSJed Brown     test:
595c4762a1bSJed Brown       requires: !single
596cbb17d71SDavid Wells       suffix: fischer_guess_3
597cbb17d71SDavid Wells       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 3,10 -pc_type none -ksp_converged_reason
598cbb17d71SDavid Wells 
599cbb17d71SDavid Wells     test:
600cbb17d71SDavid Wells       requires: !single
601c4762a1bSJed Brown       suffix: stringview
602c4762a1bSJed Brown       args: -nox -ts_type rosw -test_string_viewer
603c4762a1bSJed Brown 
604c4762a1bSJed Brown     test:
605c4762a1bSJed Brown       requires: !single
606c4762a1bSJed Brown       suffix: stringview_euler
607c4762a1bSJed Brown       args: -nox -ts_type euler -test_string_viewer
608c4762a1bSJed Brown 
609c4762a1bSJed Brown TEST*/
610