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