Actual source code: ex6.c
2: /* Program usage: ex3 [-help] [all PETSc options] */
4: static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\
5: Input parameters include:\n\
6: -m <points>, where <points> = number of grid points\n\
7: -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
8: -debug : Activate debugging printouts\n\
9: -nox : Deactivate x-window graphics\n\n";
11: /*
12: Concepts: TS^time-dependent linear problems
13: Concepts: TS^heat equation
14: Concepts: TS^diffusion equation
15: Routines: TSCreate(); TSSetSolution(); TSSetRHSMatrix();
16: Routines: TSSetInitialTimeStep(); TSSetDuration(); TSMonitorSet();
17: Routines: TSSetFromOptions(); TSStep(); TSDestroy();
18: Routines: TSSetTimeStep(); TSGetTimeStep();
19: Processors: 1
20: */
22: /* ------------------------------------------------------------------------
24: This program solves the one-dimensional heat equation (also called the
25: diffusion equation),
26: u_t = u_xx,
27: on the domain 0 <= x <= 1, with the boundary conditions
28: u(t,0) = 0, u(t,1) = 0,
29: and the initial condition
30: u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
31: This is a linear, second-order, parabolic equation.
33: We discretize the right-hand side using finite differences with
34: uniform grid spacing h:
35: u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
36: We then demonstrate time evolution using the various TS methods by
37: running the program via
38: ex3 -ts_type <timestepping solver>
40: We compare the approximate solution with the exact solution, given by
41: u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
42: 3*exp(-4*pi*pi*t) * sin(2*pi*x)
44: Notes:
45: This code demonstrates the TS solver interface to two variants of
46: linear problems, u_t = f(u,t), namely
47: - time-dependent f: f(u,t) is a function of t
48: - time-independent f: f(u,t) is simply f(u)
50: The parallel version of this code is ts/examples/tutorials/ex4.c
52: ------------------------------------------------------------------------- */
54: /*
55: Include "ts.h" so that we can use TS solvers. Note that this file
56: automatically includes:
57: petsc.h - base PETSc routines vec.h - vectors
58: sys.h - system routines mat.h - matrices
59: is.h - index sets ksp.h - Krylov subspace methods
60: viewer.h - viewers pc.h - preconditioners
61: snes.h - nonlinear solvers
62: */
64: #include petscts.h
66: /*
67: User-defined application context - contains data needed by the
68: application-provided call-back routines.
69: */
70: typedef struct {
71: Vec solution; /* global exact solution vector */
72: PetscInt m; /* total number of grid points */
73: PetscReal h; /* mesh width h = 1/(m-1) */
74: PetscTruth debug; /* flag (1 indicates activation of debugging printouts) */
75: PetscViewer viewer1, viewer2; /* viewers for the solution and error */
76: PetscReal norm_2, norm_max; /* error norms */
77: } AppCtx;
79: /*
80: User-defined routines
81: */
90: int main(int argc,char **argv)
91: {
92: AppCtx appctx; /* user-defined application context */
93: TS ts; /* timestepping context */
94: Mat A; /* matrix data structure */
95: Vec u; /* approximate solution vector */
96: PetscReal time_total_max = 100.0; /* default max total time */
97: PetscInt time_steps_max = 100; /* default max timesteps */
98: PetscDraw draw; /* drawing context */
100: PetscInt steps, m;
101: PetscMPIInt size;
102: PetscReal dt;
103: PetscReal ftime;
104: PetscTruth flg;
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Initialize program and set problem parameters
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108:
109: PetscInitialize(&argc,&argv,(char*)0,help);
110: MPI_Comm_size(PETSC_COMM_WORLD,&size);
111: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
113: m = 60;
114: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
115: PetscOptionsHasName(PETSC_NULL,"-debug",&appctx.debug);
116: appctx.m = m;
117: appctx.h = 1.0/(m-1.0);
118: appctx.norm_2 = 0.0;
119: appctx.norm_max = 0.0;
120: PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");
122: PetscOptionsGetInt(PETSC_NULL,"-time_steps_max",&time_steps_max,PETSC_NULL);
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Create vector data structures
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: /*
129: Create vector data structures for approximate and exact solutions
130: */
131: VecCreateSeq(PETSC_COMM_SELF,m,&u);
132: VecDuplicate(u,&appctx.solution);
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Set up displays to show graphs of the solution and error
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
139: PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
140: PetscDrawSetDoubleBuffer(draw);
141: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
142: PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
143: PetscDrawSetDoubleBuffer(draw);
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Create timestepping solver context
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: TSCreate(PETSC_COMM_SELF,&ts);
150: TSSetProblemType(ts,TS_LINEAR);
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Set optional user-defined monitoring routine
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: TSMonitorSet(ts,Monitor,&appctx,PETSC_NULL);
158: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: Create matrix data structure; set matrix evaluation routine.
161: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: MatCreate(PETSC_COMM_SELF,&A);
164: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
165: MatSetFromOptions(A);
167: PetscOptionsHasName(PETSC_NULL,"-time_dependent_rhs",&flg);
168: if (flg) {
169: /*
170: For linear problems with a time-dependent f(u,t) in the equation
171: u_t = f(u,t), the user provides the discretized right-hand-side
172: as a time-dependent matrix.
173: */
174: TSSetRHSMatrix(ts,A,A,RHSMatrixHeat,&appctx);
175: } else {
176: /*
177: For linear problems with a time-independent f(u) in the equation
178: u_t = f(u), the user provides the discretized right-hand-side
179: as a matrix only once, and then sets a null matrix evaluation
180: routine.
181: */
182: MatStructure A_structure;
183: RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);
184: TSSetRHSMatrix(ts,A,A,PETSC_NULL,&appctx);
185: }
187: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188: Set solution vector and initial timestep
189: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191: dt = appctx.h*appctx.h/2.0;
192: TSSetInitialTimeStep(ts,0.0,dt);
193: TSSetSolution(ts,u);
195: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196: Customize timestepping solver:
197: - Set the solution method to be the Backward Euler method.
198: - Set timestepping duration info
199: Then set runtime options, which can override these defaults.
200: For example,
201: -ts_max_steps <maxsteps> -ts_max_time <maxtime>
202: to override the defaults set by TSSetDuration().
203: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205: TSSetDuration(ts,time_steps_max,time_total_max);
206: TSSetFromOptions(ts);
208: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209: Solve the problem
210: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212: /*
213: Evaluate initial conditions
214: */
215: InitialConditions(u,&appctx);
217: /*
218: Run the timestepping solver
219: */
220: TSStep(ts,&steps,&ftime);
222: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223: View timestepping solver info
224: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226: PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %G, avg. error (max norm) = %G\n",
227: appctx.norm_2/steps,appctx.norm_max/steps);
228: TSView(ts,PETSC_VIEWER_STDOUT_SELF);
230: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231: Free work space. All PETSc objects should be destroyed when they
232: are no longer needed.
233: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
235: TSDestroy(ts);
236: MatDestroy(A);
237: VecDestroy(u);
238: PetscViewerDestroy(appctx.viewer1);
239: PetscViewerDestroy(appctx.viewer2);
240: VecDestroy(appctx.solution);
242: /*
243: Always call PetscFinalize() before exiting a program. This routine
244: - finalizes the PETSc libraries as well as MPI
245: - provides summary and diagnostic information if certain runtime
246: options are chosen (e.g., -log_summary).
247: */
248: PetscFinalize();
249: return 0;
250: }
251: /* --------------------------------------------------------------------- */
254: /*
255: InitialConditions - Computes the solution at the initial time.
257: Input Parameter:
258: u - uninitialized solution vector (global)
259: appctx - user-defined application context
261: Output Parameter:
262: u - vector with solution at initial time (global)
263: */
264: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
265: {
266: PetscScalar *u_localptr;
267: PetscInt i;
270: /*
271: Get a pointer to vector data.
272: - For default PETSc vectors, VecGetArray() returns a pointer to
273: the data array. Otherwise, the routine is implementation dependent.
274: - You MUST call VecRestoreArray() when you no longer need access to
275: the array.
276: - Note that the Fortran interface to VecGetArray() differs from the
277: C version. See the users manual for details.
278: */
279: VecGetArray(u,&u_localptr);
281: /*
282: We initialize the solution array by simply writing the solution
283: directly into the array locations. Alternatively, we could use
284: VecSetValues() or VecSetValuesLocal().
285: */
286: for (i=0; i<appctx->m; i++) {
287: u_localptr[i] = sin(PETSC_PI*i*6.*appctx->h) + 3.*sin(PETSC_PI*i*2.*appctx->h);
288: }
290: /*
291: Restore vector
292: */
293: VecRestoreArray(u,&u_localptr);
295: /*
296: Print debugging information if desired
297: */
298: if (appctx->debug) {
299: printf("initial guess vector\n");
300: VecView(u,PETSC_VIEWER_STDOUT_SELF);
301: }
303: return 0;
304: }
305: /* --------------------------------------------------------------------- */
308: /*
309: ExactSolution - Computes the exact solution at a given time.
311: Input Parameters:
312: t - current time
313: solution - vector in which exact solution will be computed
314: appctx - user-defined application context
316: Output Parameter:
317: solution - vector with the newly computed exact solution
318: */
319: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
320: {
321: PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
322: PetscInt i;
325: /*
326: Get a pointer to vector data.
327: */
328: VecGetArray(solution,&s_localptr);
330: /*
331: Simply write the solution directly into the array locations.
332: Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
333: */
334: ex1 = exp(-36.*PETSC_PI*PETSC_PI*t); ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
335: sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h;
336: for (i=0; i<appctx->m; i++) {
337: s_localptr[i] = sin(PetscRealPart(sc1)*(PetscReal)i)*ex1 + 3.*sin(PetscRealPart(sc2)*(PetscReal)i)*ex2;
338: }
340: /*
341: Restore vector
342: */
343: VecRestoreArray(solution,&s_localptr);
344: return 0;
345: }
346: /* --------------------------------------------------------------------- */
349: /*
350: Monitor - User-provided routine to monitor the solution computed at
351: each timestep. This example plots the solution and computes the
352: error in two different norms.
354: This example also demonstrates changing the timestep via TSSetTimeStep().
356: Input Parameters:
357: ts - the timestep context
358: step - the count of the current step (with 0 meaning the
359: initial condition)
360: crtime - the current time
361: u - the solution at this timestep
362: ctx - the user-provided context for this monitoring routine.
363: In this case we use the application context which contains
364: information about the problem size, workspace and the exact
365: solution.
366: */
367: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal crtime,Vec u,void *ctx)
368: {
369: AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */
371: PetscReal norm_2, norm_max, dt, dttol;
372: PetscTruth flg;
374: /*
375: View a graph of the current iterate
376: */
377: VecView(u,appctx->viewer2);
379: /*
380: Compute the exact solution
381: */
382: ExactSolution(crtime,appctx->solution,appctx);
384: /*
385: Print debugging information if desired
386: */
387: if (appctx->debug) {
388: printf("Computed solution vector\n");
389: VecView(u,PETSC_VIEWER_STDOUT_SELF);
390: printf("Exact solution vector\n");
391: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
392: }
394: /*
395: Compute the 2-norm and max-norm of the error
396: */
397: VecAXPY(appctx->solution,-1.0,u);
398: VecNorm(appctx->solution,NORM_2,&norm_2);
399: norm_2 = sqrt(appctx->h)*norm_2;
400: VecNorm(appctx->solution,NORM_MAX,&norm_max);
402: TSGetTimeStep(ts,&dt);
403: printf("Timestep %d: step size = %G, time = %G, 2-norm error = %G, max norm error = %G\n",
404: (int)step,dt,crtime,norm_2,norm_max);
405: appctx->norm_2 += norm_2;
406: appctx->norm_max += norm_max;
408: dttol = .0001;
409: PetscOptionsGetReal(PETSC_NULL,"-dttol",&dttol,&flg);
410: if (dt < dttol) {
411: dt *= .999;
412: TSSetTimeStep(ts,dt);
413: }
415: /*
416: View a graph of the error
417: */
418: VecView(appctx->solution,appctx->viewer1);
420: /*
421: Print debugging information if desired
422: */
423: if (appctx->debug) {
424: printf("Error vector\n");
425: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
426: }
428: return 0;
429: }
430: /* --------------------------------------------------------------------- */
433: /*
434: RHSMatrixHeat - User-provided routine to compute the right-hand-side
435: matrix for the heat equation.
437: Input Parameters:
438: ts - the TS context
439: t - current time
440: global_in - global input vector
441: dummy - optional user-defined context, as set by TSetRHSJacobian()
443: Output Parameters:
444: AA - Jacobian matrix
445: BB - optionally different preconditioning matrix
446: str - flag indicating matrix structure
448: Notes:
449: Recall that MatSetValues() uses 0-based row and column numbers
450: in Fortran as well as in C.
451: */
452: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
453: {
454: Mat A = *AA; /* Jacobian matrix */
455: AppCtx *appctx = (AppCtx *) ctx; /* user-defined application context */
456: PetscInt mstart = 0;
457: PetscInt mend = appctx->m;
459: PetscInt i, idx[3];
460: PetscScalar v[3], stwo = -2./(appctx->h*appctx->h), sone = -.5*stwo;
462: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
463: Compute entries for the locally owned part of the matrix
464: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
465: /*
466: Set matrix rows corresponding to boundary data
467: */
469: mstart = 0;
470: v[0] = 1.0;
471: MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
472: mstart++;
474: mend--;
475: v[0] = 1.0;
476: MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
478: /*
479: Set matrix rows corresponding to interior data. We construct the
480: matrix one row at a time.
481: */
482: v[0] = sone; v[1] = stwo; v[2] = sone;
483: for ( i=mstart; i<mend; i++ ) {
484: idx[0] = i-1; idx[1] = i; idx[2] = i+1;
485: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
486: }
488: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
489: Complete the matrix assembly process and set some options
490: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
491: /*
492: Assemble matrix, using the 2-step process:
493: MatAssemblyBegin(), MatAssemblyEnd()
494: Computations can be done while messages are in transition
495: by placing code between these two statements.
496: */
497: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
498: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
500: /*
501: Set flag to indicate that the Jacobian matrix retains an identical
502: nonzero structure throughout all timestepping iterations (although the
503: values of the entries change). Thus, we can save some work in setting
504: up the preconditioner (e.g., no need to redo symbolic factorization for
505: ILU/ICC preconditioners).
506: - If the nonzero structure of the matrix is different during
507: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
508: must be used instead. If you are unsure whether the matrix
509: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
510: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
511: believes your assertion and does not check the structure
512: of the matrix. If you erroneously claim that the structure
513: is the same when it actually is not, the new preconditioner
514: will not function correctly. Thus, use this optimization
515: feature with caution!
516: */
517: *str = SAME_NONZERO_PATTERN;
519: /*
520: Set and option to indicate that we will never add a new nonzero location
521: to the matrix. If we do, it will generate an error.
522: */
523: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);
525: return 0;
526: }
527: /* --------------------------------------------------------------------- */
530: /*
531: Input Parameters:
532: ts - the TS context
533: t - current time
534: f - function
535: ctx - optional user-defined context, as set by TSetBCFunction()
536: */
537: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
538: {
539: AppCtx *appctx = (AppCtx *) ctx; /* user-defined application context */
541: PetscInt m = appctx->m;
542: PetscScalar *fa;
544: VecGetArray(f,&fa);
545: fa[0] = 0.0;
546: fa[m-1] = 1.0;
547: VecRestoreArray(f,&fa);
548: printf("t=%g\n",t);
549:
550: return 0;
551: }