xref: /petsc/src/ts/tutorials/ex4.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: n
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) = 0, u(t,1) = 0,
23*c4762a1bSJed Brown    and the initial condition
24*c4762a1bSJed Brown        u(0,x) = sin(6*pi*x) + 3*sin(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        mpiexec -n <procs> 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) * sin(6*pi*x) +
36*c4762a1bSJed Brown                       3*exp(-4*pi*pi*t) * sin(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 f(u)
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown     The uniprocessor version of this code is ts/tutorials/ex3.c
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown   ------------------------------------------------------------------------- */
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown /*
49*c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs) to manage
50*c4762a1bSJed Brown    the parallel grid.  Include "petscts.h" so that we can use TS solvers.
51*c4762a1bSJed Brown    Note that this file automatically includes:
52*c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h  - vectors
53*c4762a1bSJed Brown      petscmat.h  - matrices
54*c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
55*c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h   - preconditioners
56*c4762a1bSJed Brown      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
57*c4762a1bSJed Brown */
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown #include <petscdm.h>
60*c4762a1bSJed Brown #include <petscdmda.h>
61*c4762a1bSJed Brown #include <petscts.h>
62*c4762a1bSJed Brown #include <petscdraw.h>
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown /*
65*c4762a1bSJed Brown    User-defined application context - contains data needed by the
66*c4762a1bSJed Brown    application-provided call-back routines.
67*c4762a1bSJed Brown */
68*c4762a1bSJed Brown typedef struct {
69*c4762a1bSJed Brown   MPI_Comm    comm;              /* communicator */
70*c4762a1bSJed Brown   DM          da;                /* distributed array data structure */
71*c4762a1bSJed Brown   Vec         localwork;         /* local ghosted work vector */
72*c4762a1bSJed Brown   Vec         u_local;           /* local ghosted approximate solution vector */
73*c4762a1bSJed Brown   Vec         solution;          /* global exact solution vector */
74*c4762a1bSJed Brown   PetscInt    m;                 /* total number of grid points */
75*c4762a1bSJed Brown   PetscReal   h;                 /* mesh width h = 1/(m-1) */
76*c4762a1bSJed Brown   PetscBool   debug;             /* flag (1 indicates activation of debugging printouts) */
77*c4762a1bSJed Brown   PetscViewer viewer1,viewer2;  /* viewers for the solution and error */
78*c4762a1bSJed Brown   PetscReal   norm_2,norm_max;  /* error norms */
79*c4762a1bSJed Brown } AppCtx;
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown /*
82*c4762a1bSJed Brown    User-defined routines
83*c4762a1bSJed Brown */
84*c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
85*c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
86*c4762a1bSJed Brown extern PetscErrorCode RHSFunctionHeat(TS,PetscReal,Vec,Vec,void*);
87*c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
88*c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown int main(int argc,char **argv)
91*c4762a1bSJed Brown {
92*c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
93*c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
94*c4762a1bSJed Brown   Mat            A;                      /* matrix data structure */
95*c4762a1bSJed Brown   Vec            u;                      /* approximate solution vector */
96*c4762a1bSJed Brown   PetscReal      time_total_max = 1.0;   /* default max total time */
97*c4762a1bSJed Brown   PetscInt       time_steps_max = 100;   /* default max timesteps */
98*c4762a1bSJed Brown   PetscDraw      draw;                   /* drawing context */
99*c4762a1bSJed Brown   PetscErrorCode ierr;
100*c4762a1bSJed Brown   PetscInt       steps,m;
101*c4762a1bSJed Brown   PetscMPIInt    size;
102*c4762a1bSJed Brown   PetscReal      dt,ftime;
103*c4762a1bSJed Brown   PetscBool      flg;
104*c4762a1bSJed Brown   TSProblemType  tsproblem = TS_LINEAR;
105*c4762a1bSJed Brown 
106*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107*c4762a1bSJed Brown      Initialize program and set problem parameters
108*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   ierr        = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
111*c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown   m               = 60;
114*c4762a1bSJed Brown   ierr            = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr            = PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);CHKERRQ(ierr);
116*c4762a1bSJed Brown   appctx.m        = m;
117*c4762a1bSJed Brown   appctx.h        = 1.0/(m-1.0);
118*c4762a1bSJed Brown   appctx.norm_2   = 0.0;
119*c4762a1bSJed Brown   appctx.norm_max = 0.0;
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Solving a linear TS problem, number of processors = %d\n",size);CHKERRQ(ierr);
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125*c4762a1bSJed Brown      Create vector data structures
126*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127*c4762a1bSJed Brown   /*
128*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
129*c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
130*c4762a1bSJed Brown      total grid values spread equally among all the processors.
131*c4762a1bSJed Brown   */
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,m,1,1,NULL,&appctx.da);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = DMSetFromOptions(appctx.da);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = DMSetUp(appctx.da);CHKERRQ(ierr);
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown   /*
138*c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
139*c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
140*c4762a1bSJed Brown      have the same types.
141*c4762a1bSJed Brown   */
142*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(appctx.da,&u);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = DMCreateLocalVector(appctx.da,&appctx.u_local);CHKERRQ(ierr);
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown   /*
146*c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
147*c4762a1bSJed Brown      create global work vector for storing exact solution.
148*c4762a1bSJed Brown   */
149*c4762a1bSJed Brown   ierr = VecDuplicate(appctx.u_local,&appctx.localwork);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = VecDuplicate(u,&appctx.solution);CHKERRQ(ierr);
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153*c4762a1bSJed Brown      Set up displays to show graphs of the solution and error
154*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,380,400,160,&appctx.viewer1);CHKERRQ(ierr);
157*c4762a1bSJed Brown   ierr = PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,0,400,160,&appctx.viewer2);CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr);
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164*c4762a1bSJed Brown      Create timestepping solver context
165*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
168*c4762a1bSJed Brown 
169*c4762a1bSJed Brown   flg  = PETSC_FALSE;
170*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-nonlinear",&flg,NULL);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,flg ? TS_NONLINEAR : TS_LINEAR);CHKERRQ(ierr);
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174*c4762a1bSJed Brown      Set optional user-defined monitoring routine
175*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176*c4762a1bSJed Brown   ierr = TSMonitorSet(ts,Monitor,&appctx,NULL);CHKERRQ(ierr);
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown      Create matrix data structure; set matrix evaluation routine.
181*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   flg  = PETSC_FALSE;
189*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL);CHKERRQ(ierr);
190*c4762a1bSJed Brown   if (flg) {
191*c4762a1bSJed Brown     /*
192*c4762a1bSJed Brown        For linear problems with a time-dependent f(u,t) in the equation
193*c4762a1bSJed Brown        u_t = f(u,t), the user provides the discretized right-hand-side
194*c4762a1bSJed Brown        as a time-dependent matrix.
195*c4762a1bSJed Brown     */
196*c4762a1bSJed Brown     ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr);
197*c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);CHKERRQ(ierr);
198*c4762a1bSJed Brown   } else {
199*c4762a1bSJed Brown     /*
200*c4762a1bSJed Brown        For linear problems with a time-independent f(u) in the equation
201*c4762a1bSJed Brown        u_t = f(u), the user provides the discretized right-hand-side
202*c4762a1bSJed Brown        as a matrix only once, and then sets a null matrix evaluation
203*c4762a1bSJed Brown        routine.
204*c4762a1bSJed Brown     */
205*c4762a1bSJed Brown     ierr = RHSMatrixHeat(ts,0.0,u,A,A,&appctx);CHKERRQ(ierr);
206*c4762a1bSJed Brown     ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr);
207*c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);CHKERRQ(ierr);
208*c4762a1bSJed Brown   }
209*c4762a1bSJed Brown 
210*c4762a1bSJed Brown   if (tsproblem == TS_NONLINEAR) {
211*c4762a1bSJed Brown     SNES snes;
212*c4762a1bSJed Brown     ierr = TSSetRHSFunction(ts,NULL,RHSFunctionHeat,&appctx);CHKERRQ(ierr);
213*c4762a1bSJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
214*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefault,NULL);CHKERRQ(ierr);
215*c4762a1bSJed Brown   }
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218*c4762a1bSJed Brown      Set solution vector and initial timestep
219*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown   dt   = appctx.h*appctx.h/2.0;
222*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
223*c4762a1bSJed Brown   ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
224*c4762a1bSJed Brown 
225*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226*c4762a1bSJed Brown      Customize timestepping solver:
227*c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
228*c4762a1bSJed Brown        - Set timestepping duration info
229*c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
230*c4762a1bSJed Brown      For example,
231*c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
232*c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
233*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
234*c4762a1bSJed Brown 
235*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,time_steps_max);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,time_total_max);CHKERRQ(ierr);
237*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
239*c4762a1bSJed Brown 
240*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241*c4762a1bSJed Brown      Solve the problem
242*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243*c4762a1bSJed Brown 
244*c4762a1bSJed Brown   /*
245*c4762a1bSJed Brown      Evaluate initial conditions
246*c4762a1bSJed Brown   */
247*c4762a1bSJed Brown   ierr = InitialConditions(u,&appctx);CHKERRQ(ierr);
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown   /*
250*c4762a1bSJed Brown      Run the timestepping solver
251*c4762a1bSJed Brown   */
252*c4762a1bSJed Brown   ierr = TSSolve(ts,u);CHKERRQ(ierr);
253*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
254*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257*c4762a1bSJed Brown      View timestepping solver info
258*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
259*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Total timesteps %D, Final time %g\n",steps,(double)ftime);CHKERRQ(ierr);
260*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Avg. error (2 norm) = %g Avg. error (max norm) = %g\n",(double)(appctx.norm_2/steps),(double)(appctx.norm_max/steps));CHKERRQ(ierr);
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
264*c4762a1bSJed Brown      are no longer needed.
265*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
266*c4762a1bSJed Brown 
267*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
270*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&appctx.viewer1);CHKERRQ(ierr);
271*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&appctx.viewer2);CHKERRQ(ierr);
272*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.localwork);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.solution);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.u_local);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = DMDestroy(&appctx.da);CHKERRQ(ierr);
276*c4762a1bSJed Brown 
277*c4762a1bSJed Brown   /*
278*c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
279*c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
280*c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
281*c4762a1bSJed Brown          options are chosen (e.g., -log_view).
282*c4762a1bSJed Brown   */
283*c4762a1bSJed Brown   ierr = PetscFinalize();
284*c4762a1bSJed Brown   return ierr;
285*c4762a1bSJed Brown }
286*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
287*c4762a1bSJed Brown /*
288*c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
289*c4762a1bSJed Brown 
290*c4762a1bSJed Brown    Input Parameter:
291*c4762a1bSJed Brown    u - uninitialized solution vector (global)
292*c4762a1bSJed Brown    appctx - user-defined application context
293*c4762a1bSJed Brown 
294*c4762a1bSJed Brown    Output Parameter:
295*c4762a1bSJed Brown    u - vector with solution at initial time (global)
296*c4762a1bSJed Brown */
297*c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
298*c4762a1bSJed Brown {
299*c4762a1bSJed Brown   PetscScalar    *u_localptr,h = appctx->h;
300*c4762a1bSJed Brown   PetscInt       i,mybase,myend;
301*c4762a1bSJed Brown   PetscErrorCode ierr;
302*c4762a1bSJed Brown 
303*c4762a1bSJed Brown   /*
304*c4762a1bSJed Brown      Determine starting point of each processor's range of
305*c4762a1bSJed Brown      grid values.
306*c4762a1bSJed Brown   */
307*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(u,&mybase,&myend);CHKERRQ(ierr);
308*c4762a1bSJed Brown 
309*c4762a1bSJed Brown   /*
310*c4762a1bSJed Brown     Get a pointer to vector data.
311*c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
312*c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
313*c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
314*c4762a1bSJed Brown       the array.
315*c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
316*c4762a1bSJed Brown       C version.  See the users manual for details.
317*c4762a1bSJed Brown   */
318*c4762a1bSJed Brown   ierr = VecGetArray(u,&u_localptr);CHKERRQ(ierr);
319*c4762a1bSJed Brown 
320*c4762a1bSJed Brown   /*
321*c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
322*c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
323*c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
324*c4762a1bSJed Brown   */
325*c4762a1bSJed Brown   for (i=mybase; i<myend; i++) u_localptr[i-mybase] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
326*c4762a1bSJed Brown 
327*c4762a1bSJed Brown   /*
328*c4762a1bSJed Brown      Restore vector
329*c4762a1bSJed Brown   */
330*c4762a1bSJed Brown   ierr = VecRestoreArray(u,&u_localptr);CHKERRQ(ierr);
331*c4762a1bSJed Brown 
332*c4762a1bSJed Brown   /*
333*c4762a1bSJed Brown      Print debugging information if desired
334*c4762a1bSJed Brown   */
335*c4762a1bSJed Brown   if (appctx->debug) {
336*c4762a1bSJed Brown     ierr = PetscPrintf(appctx->comm,"initial guess vector\n");CHKERRQ(ierr);
337*c4762a1bSJed Brown     ierr = VecView(u,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
338*c4762a1bSJed Brown   }
339*c4762a1bSJed Brown 
340*c4762a1bSJed Brown   return 0;
341*c4762a1bSJed Brown }
342*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
343*c4762a1bSJed Brown /*
344*c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
345*c4762a1bSJed Brown 
346*c4762a1bSJed Brown    Input Parameters:
347*c4762a1bSJed Brown    t - current time
348*c4762a1bSJed Brown    solution - vector in which exact solution will be computed
349*c4762a1bSJed Brown    appctx - user-defined application context
350*c4762a1bSJed Brown 
351*c4762a1bSJed Brown    Output Parameter:
352*c4762a1bSJed Brown    solution - vector with the newly computed exact solution
353*c4762a1bSJed Brown */
354*c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
355*c4762a1bSJed Brown {
356*c4762a1bSJed Brown   PetscScalar    *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2;
357*c4762a1bSJed Brown   PetscInt       i,mybase,myend;
358*c4762a1bSJed Brown   PetscErrorCode ierr;
359*c4762a1bSJed Brown 
360*c4762a1bSJed Brown   /*
361*c4762a1bSJed Brown      Determine starting and ending points of each processor's
362*c4762a1bSJed Brown      range of grid values
363*c4762a1bSJed Brown   */
364*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(solution,&mybase,&myend);CHKERRQ(ierr);
365*c4762a1bSJed Brown 
366*c4762a1bSJed Brown   /*
367*c4762a1bSJed Brown      Get a pointer to vector data.
368*c4762a1bSJed Brown   */
369*c4762a1bSJed Brown   ierr = VecGetArray(solution,&s_localptr);CHKERRQ(ierr);
370*c4762a1bSJed Brown 
371*c4762a1bSJed Brown   /*
372*c4762a1bSJed Brown      Simply write the solution directly into the array locations.
373*c4762a1bSJed Brown      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
374*c4762a1bSJed Brown   */
375*c4762a1bSJed Brown   ex1 = PetscExpReal(-36.*PETSC_PI*PETSC_PI*t); ex2 = PetscExpReal(-4.*PETSC_PI*PETSC_PI*t);
376*c4762a1bSJed Brown   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
377*c4762a1bSJed Brown   for (i=mybase; i<myend; i++) s_localptr[i-mybase] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
378*c4762a1bSJed Brown 
379*c4762a1bSJed Brown   /*
380*c4762a1bSJed Brown      Restore vector
381*c4762a1bSJed Brown   */
382*c4762a1bSJed Brown   ierr = VecRestoreArray(solution,&s_localptr);CHKERRQ(ierr);
383*c4762a1bSJed Brown   return 0;
384*c4762a1bSJed Brown }
385*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
386*c4762a1bSJed Brown /*
387*c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
388*c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
389*c4762a1bSJed Brown    error in two different norms.
390*c4762a1bSJed Brown 
391*c4762a1bSJed Brown    Input Parameters:
392*c4762a1bSJed Brown    ts     - the timestep context
393*c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
394*c4762a1bSJed Brown              initial condition)
395*c4762a1bSJed Brown    time   - the current time
396*c4762a1bSJed Brown    u      - the solution at this timestep
397*c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
398*c4762a1bSJed Brown             In this case we use the application context which contains
399*c4762a1bSJed Brown             information about the problem size, workspace and the exact
400*c4762a1bSJed Brown             solution.
401*c4762a1bSJed Brown */
402*c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
403*c4762a1bSJed Brown {
404*c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
405*c4762a1bSJed Brown   PetscErrorCode ierr;
406*c4762a1bSJed Brown   PetscReal      norm_2,norm_max;
407*c4762a1bSJed Brown 
408*c4762a1bSJed Brown   /*
409*c4762a1bSJed Brown      View a graph of the current iterate
410*c4762a1bSJed Brown   */
411*c4762a1bSJed Brown   ierr = VecView(u,appctx->viewer2);CHKERRQ(ierr);
412*c4762a1bSJed Brown 
413*c4762a1bSJed Brown   /*
414*c4762a1bSJed Brown      Compute the exact solution
415*c4762a1bSJed Brown   */
416*c4762a1bSJed Brown   ierr = ExactSolution(time,appctx->solution,appctx);CHKERRQ(ierr);
417*c4762a1bSJed Brown 
418*c4762a1bSJed Brown   /*
419*c4762a1bSJed Brown      Print debugging information if desired
420*c4762a1bSJed Brown   */
421*c4762a1bSJed Brown   if (appctx->debug) {
422*c4762a1bSJed Brown     ierr = PetscPrintf(appctx->comm,"Computed solution vector\n");CHKERRQ(ierr);
423*c4762a1bSJed Brown     ierr = VecView(u,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
424*c4762a1bSJed Brown     ierr = PetscPrintf(appctx->comm,"Exact solution vector\n");CHKERRQ(ierr);
425*c4762a1bSJed Brown     ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
426*c4762a1bSJed Brown   }
427*c4762a1bSJed Brown 
428*c4762a1bSJed Brown   /*
429*c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
430*c4762a1bSJed Brown   */
431*c4762a1bSJed Brown   ierr   = VecAXPY(appctx->solution,-1.0,u);CHKERRQ(ierr);
432*c4762a1bSJed Brown   ierr   = VecNorm(appctx->solution,NORM_2,&norm_2);CHKERRQ(ierr);
433*c4762a1bSJed Brown   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
434*c4762a1bSJed Brown   ierr   = VecNorm(appctx->solution,NORM_MAX,&norm_max);CHKERRQ(ierr);
435*c4762a1bSJed Brown   if (norm_2   < 1e-14) norm_2   = 0;
436*c4762a1bSJed Brown   if (norm_max < 1e-14) norm_max = 0;
437*c4762a1bSJed Brown 
438*c4762a1bSJed Brown   /*
439*c4762a1bSJed Brown      PetscPrintf() causes only the first processor in this
440*c4762a1bSJed Brown      communicator to print the timestep information.
441*c4762a1bSJed Brown   */
442*c4762a1bSJed Brown   ierr = PetscPrintf(appctx->comm,"Timestep %D: time = %g 2-norm error = %g max norm error = %g\n",step,(double)time,(double)norm_2,(double)norm_max);CHKERRQ(ierr);
443*c4762a1bSJed Brown   appctx->norm_2   += norm_2;
444*c4762a1bSJed Brown   appctx->norm_max += norm_max;
445*c4762a1bSJed Brown 
446*c4762a1bSJed Brown   /*
447*c4762a1bSJed Brown      View a graph of the error
448*c4762a1bSJed Brown   */
449*c4762a1bSJed Brown   ierr = VecView(appctx->solution,appctx->viewer1);CHKERRQ(ierr);
450*c4762a1bSJed Brown 
451*c4762a1bSJed Brown   /*
452*c4762a1bSJed Brown      Print debugging information if desired
453*c4762a1bSJed Brown   */
454*c4762a1bSJed Brown   if (appctx->debug) {
455*c4762a1bSJed Brown     ierr = PetscPrintf(appctx->comm,"Error vector\n");CHKERRQ(ierr);
456*c4762a1bSJed Brown     ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
457*c4762a1bSJed Brown   }
458*c4762a1bSJed Brown 
459*c4762a1bSJed Brown   return 0;
460*c4762a1bSJed Brown }
461*c4762a1bSJed Brown 
462*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
463*c4762a1bSJed Brown /*
464*c4762a1bSJed Brown    RHSMatrixHeat - User-provided routine to compute the right-hand-side
465*c4762a1bSJed Brown    matrix for the heat equation.
466*c4762a1bSJed Brown 
467*c4762a1bSJed Brown    Input Parameters:
468*c4762a1bSJed Brown    ts - the TS context
469*c4762a1bSJed Brown    t - current time
470*c4762a1bSJed Brown    global_in - global input vector
471*c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
472*c4762a1bSJed Brown 
473*c4762a1bSJed Brown    Output Parameters:
474*c4762a1bSJed Brown    AA - Jacobian matrix
475*c4762a1bSJed Brown    BB - optionally different preconditioning matrix
476*c4762a1bSJed Brown    str - flag indicating matrix structure
477*c4762a1bSJed Brown 
478*c4762a1bSJed Brown   Notes:
479*c4762a1bSJed Brown   RHSMatrixHeat computes entries for the locally owned part of the system.
480*c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
481*c4762a1bSJed Brown      contiguous chunks of rows across the processors.
482*c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
483*c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
484*c4762a1bSJed Brown      appropriate processor during matrix assembly).
485*c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
486*c4762a1bSJed Brown      using MatSetValues(); we could alternatively use MatSetValuesLocal().
487*c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
488*c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
489*c4762a1bSJed Brown      in Fortran as well as in C.
490*c4762a1bSJed Brown */
491*c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
492*c4762a1bSJed Brown {
493*c4762a1bSJed Brown   Mat            A       = AA;              /* Jacobian matrix */
494*c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
495*c4762a1bSJed Brown   PetscErrorCode ierr;
496*c4762a1bSJed Brown   PetscInt       i,mstart,mend,idx[3];
497*c4762a1bSJed Brown   PetscScalar    v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;
498*c4762a1bSJed Brown 
499*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
500*c4762a1bSJed Brown      Compute entries for the locally owned part of the matrix
501*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
502*c4762a1bSJed Brown 
503*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&mstart,&mend);CHKERRQ(ierr);
504*c4762a1bSJed Brown 
505*c4762a1bSJed Brown   /*
506*c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
507*c4762a1bSJed Brown   */
508*c4762a1bSJed Brown 
509*c4762a1bSJed Brown   if (mstart == 0) {  /* first processor only */
510*c4762a1bSJed Brown     v[0] = 1.0;
511*c4762a1bSJed Brown     ierr = MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);CHKERRQ(ierr);
512*c4762a1bSJed Brown     mstart++;
513*c4762a1bSJed Brown   }
514*c4762a1bSJed Brown 
515*c4762a1bSJed Brown   if (mend == appctx->m) { /* last processor only */
516*c4762a1bSJed Brown     mend--;
517*c4762a1bSJed Brown     v[0] = 1.0;
518*c4762a1bSJed Brown     ierr = MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);CHKERRQ(ierr);
519*c4762a1bSJed Brown   }
520*c4762a1bSJed Brown 
521*c4762a1bSJed Brown   /*
522*c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
523*c4762a1bSJed Brown      matrix one row at a time.
524*c4762a1bSJed Brown   */
525*c4762a1bSJed Brown   v[0] = sone; v[1] = stwo; v[2] = sone;
526*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
527*c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
528*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
529*c4762a1bSJed Brown   }
530*c4762a1bSJed Brown 
531*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
532*c4762a1bSJed Brown      Complete the matrix assembly process and set some options
533*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
534*c4762a1bSJed Brown   /*
535*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
536*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
537*c4762a1bSJed Brown      Computations can be done while messages are in transition
538*c4762a1bSJed Brown      by placing code between these two statements.
539*c4762a1bSJed Brown   */
540*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
541*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
542*c4762a1bSJed Brown 
543*c4762a1bSJed Brown   /*
544*c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
545*c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
546*c4762a1bSJed Brown   */
547*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
548*c4762a1bSJed Brown 
549*c4762a1bSJed Brown   return 0;
550*c4762a1bSJed Brown }
551*c4762a1bSJed Brown 
552*c4762a1bSJed Brown PetscErrorCode RHSFunctionHeat(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
553*c4762a1bSJed Brown {
554*c4762a1bSJed Brown   PetscErrorCode ierr;
555*c4762a1bSJed Brown   Mat            A;
556*c4762a1bSJed Brown 
557*c4762a1bSJed Brown   PetscFunctionBeginUser;
558*c4762a1bSJed Brown   ierr = TSGetRHSJacobian(ts,&A,NULL,NULL,&ctx);CHKERRQ(ierr);
559*c4762a1bSJed Brown   ierr = RHSMatrixHeat(ts,t,globalin,A,NULL,ctx);CHKERRQ(ierr);
560*c4762a1bSJed Brown   /* ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
561*c4762a1bSJed Brown   ierr = MatMult(A,globalin,globalout);CHKERRQ(ierr);
562*c4762a1bSJed Brown   PetscFunctionReturn(0);
563*c4762a1bSJed Brown }
564*c4762a1bSJed Brown 
565*c4762a1bSJed Brown /*TEST
566*c4762a1bSJed Brown 
567*c4762a1bSJed Brown     test:
568*c4762a1bSJed Brown       args: -ts_view -nox
569*c4762a1bSJed Brown 
570*c4762a1bSJed Brown     test:
571*c4762a1bSJed Brown       suffix: 2
572*c4762a1bSJed Brown       args: -ts_view -nox
573*c4762a1bSJed Brown       nsize: 3
574*c4762a1bSJed Brown 
575*c4762a1bSJed Brown     test:
576*c4762a1bSJed Brown       suffix: 3
577*c4762a1bSJed Brown       args: -ts_view -nox -nonlinear
578*c4762a1bSJed Brown 
579*c4762a1bSJed Brown     test:
580*c4762a1bSJed Brown       suffix: 4
581*c4762a1bSJed Brown       args: -ts_view -nox -nonlinear
582*c4762a1bSJed Brown       nsize: 3
583*c4762a1bSJed Brown       timeoutfactor: 3
584*c4762a1bSJed Brown 
585*c4762a1bSJed Brown     test:
586*c4762a1bSJed Brown       suffix: sundials
587*c4762a1bSJed Brown       requires: sundials
588*c4762a1bSJed Brown       args: -nox -ts_type sundials -ts_max_steps 5 -nonlinear
589*c4762a1bSJed Brown       nsize: 4
590*c4762a1bSJed Brown 
591*c4762a1bSJed Brown TEST*/
592*c4762a1bSJed Brown 
593*c4762a1bSJed Brown 
594