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