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(); TSSetMatrices();
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: petscsys.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: PetscBool 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: PetscBool 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(PETSC_COMM_SELF,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: TSSetRHSFunction(ts,PETSC_NULL,TSComputeRHSFunctionLinear,&appctx);
175: TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
176: } else {
177: /*
178: For linear problems with a time-independent f(u) in the equation
179: u_t = f(u), the user provides the discretized right-hand-side
180: as a matrix only once, and then sets a null matrix evaluation
181: routine.
182: */
183: MatStructure A_structure;
184: RHSMatrixHeat(ts,0.0,u,&A,&A,&A_structure,&appctx);
185: TSSetRHSFunction(ts,PETSC_NULL,TSComputeRHSFunctionLinear,&appctx);
186: TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
187: }
189: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190: Set solution vector and initial timestep
191: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193: dt = appctx.h*appctx.h/2.0;
194: TSSetInitialTimeStep(ts,0.0,dt);
195: TSSetSolution(ts,u);
197: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: Customize timestepping solver:
199: - Set the solution method to be the Backward Euler method.
200: - Set timestepping duration info
201: Then set runtime options, which can override these defaults.
202: For example,
203: -ts_max_steps <maxsteps> -ts_max_time <maxtime>
204: to override the defaults set by TSSetDuration().
205: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207: TSSetDuration(ts,time_steps_max,time_total_max);
208: TSSetFromOptions(ts);
210: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: Solve the problem
212: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214: /*
215: Evaluate initial conditions
216: */
217: InitialConditions(u,&appctx);
219: /*
220: Run the timestepping solver
221: */
222: TSSolve(ts,u,&ftime);
223: TSGetTimeStepNumber(ts,&steps);
225: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226: View timestepping solver info
227: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229: PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %G, avg. error (max norm) = %G\n",
230: appctx.norm_2/steps,appctx.norm_max/steps);
231: TSView(ts,PETSC_VIEWER_STDOUT_SELF);
233: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234: Free work space. All PETSc objects should be destroyed when they
235: are no longer needed.
236: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238: TSDestroy(&ts);
239: MatDestroy(&A);
240: VecDestroy(&u);
241: PetscViewerDestroy(&appctx.viewer1);
242: PetscViewerDestroy(&appctx.viewer2);
243: VecDestroy(&appctx.solution);
245: /*
246: Always call PetscFinalize() before exiting a program. This routine
247: - finalizes the PETSc libraries as well as MPI
248: - provides summary and diagnostic information if certain runtime
249: options are chosen (e.g., -log_summary).
250: */
251: PetscFinalize();
252: return 0;
253: }
254: /* --------------------------------------------------------------------- */
257: /*
258: InitialConditions - Computes the solution at the initial time.
260: Input Parameter:
261: u - uninitialized solution vector (global)
262: appctx - user-defined application context
264: Output Parameter:
265: u - vector with solution at initial time (global)
266: */
267: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
268: {
269: PetscScalar *u_localptr;
270: PetscInt i;
273: /*
274: Get a pointer to vector data.
275: - For default PETSc vectors, VecGetArray() returns a pointer to
276: the data array. Otherwise, the routine is implementation dependent.
277: - You MUST call VecRestoreArray() when you no longer need access to
278: the array.
279: - Note that the Fortran interface to VecGetArray() differs from the
280: C version. See the users manual for details.
281: */
282: VecGetArray(u,&u_localptr);
284: /*
285: We initialize the solution array by simply writing the solution
286: directly into the array locations. Alternatively, we could use
287: VecSetValues() or VecSetValuesLocal().
288: */
289: for (i=0; i<appctx->m; i++) {
290: u_localptr[i] = sin(PETSC_PI*i*6.*appctx->h) + 3.*sin(PETSC_PI*i*2.*appctx->h);
291: }
293: /*
294: Restore vector
295: */
296: VecRestoreArray(u,&u_localptr);
298: /*
299: Print debugging information if desired
300: */
301: if (appctx->debug) {
302: printf("initial guess vector\n");
303: VecView(u,PETSC_VIEWER_STDOUT_SELF);
304: }
306: return 0;
307: }
308: /* --------------------------------------------------------------------- */
311: /*
312: ExactSolution - Computes the exact solution at a given time.
314: Input Parameters:
315: t - current time
316: solution - vector in which exact solution will be computed
317: appctx - user-defined application context
319: Output Parameter:
320: solution - vector with the newly computed exact solution
321: */
322: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
323: {
324: PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
325: PetscInt i;
328: /*
329: Get a pointer to vector data.
330: */
331: VecGetArray(solution,&s_localptr);
333: /*
334: Simply write the solution directly into the array locations.
335: Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
336: */
337: ex1 = exp(-36.*PETSC_PI*PETSC_PI*t); ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
338: sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h;
339: for (i=0; i<appctx->m; i++) {
340: s_localptr[i] = sin(PetscRealPart(sc1)*(PetscReal)i)*ex1 + 3.*sin(PetscRealPart(sc2)*(PetscReal)i)*ex2;
341: }
343: /*
344: Restore vector
345: */
346: VecRestoreArray(solution,&s_localptr);
347: return 0;
348: }
349: /* --------------------------------------------------------------------- */
352: /*
353: Monitor - User-provided routine to monitor the solution computed at
354: each timestep. This example plots the solution and computes the
355: error in two different norms.
357: This example also demonstrates changing the timestep via TSSetTimeStep().
359: Input Parameters:
360: ts - the timestep context
361: step - the count of the current step (with 0 meaning the
362: initial condition)
363: crtime - the current time
364: u - the solution at this timestep
365: ctx - the user-provided context for this monitoring routine.
366: In this case we use the application context which contains
367: information about the problem size, workspace and the exact
368: solution.
369: */
370: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal crtime,Vec u,void *ctx)
371: {
372: AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */
374: PetscReal norm_2, norm_max, dt, dttol;
375: PetscBool flg;
377: /*
378: View a graph of the current iterate
379: */
380: VecView(u,appctx->viewer2);
382: /*
383: Compute the exact solution
384: */
385: ExactSolution(crtime,appctx->solution,appctx);
387: /*
388: Print debugging information if desired
389: */
390: if (appctx->debug) {
391: printf("Computed solution vector\n");
392: VecView(u,PETSC_VIEWER_STDOUT_SELF);
393: printf("Exact solution vector\n");
394: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
395: }
397: /*
398: Compute the 2-norm and max-norm of the error
399: */
400: VecAXPY(appctx->solution,-1.0,u);
401: VecNorm(appctx->solution,NORM_2,&norm_2);
402: norm_2 = sqrt(appctx->h)*norm_2;
403: VecNorm(appctx->solution,NORM_MAX,&norm_max);
405: TSGetTimeStep(ts,&dt);
406: if (norm_2 > 1.e-2){
407: printf("Timestep %d: step size = %G, time = %G, 2-norm error = %G, max norm error = %G\n",
408: (int)step,dt,crtime,norm_2,norm_max);
409: }
410: appctx->norm_2 += norm_2;
411: appctx->norm_max += norm_max;
413: dttol = .0001;
414: PetscOptionsGetReal(PETSC_NULL,"-dttol",&dttol,&flg);
415: if (dt < dttol) {
416: dt *= .999;
417: TSSetTimeStep(ts,dt);
418: }
420: /*
421: View a graph of the error
422: */
423: VecView(appctx->solution,appctx->viewer1);
425: /*
426: Print debugging information if desired
427: */
428: if (appctx->debug) {
429: printf("Error vector\n");
430: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
431: }
433: return 0;
434: }
435: /* --------------------------------------------------------------------- */
438: /*
439: RHSMatrixHeat - User-provided routine to compute the right-hand-side
440: matrix for the heat equation.
442: Input Parameters:
443: ts - the TS context
444: t - current time
445: global_in - global input vector
446: dummy - optional user-defined context, as set by TSetRHSJacobian()
448: Output Parameters:
449: AA - Jacobian matrix
450: BB - optionally different preconditioning matrix
451: str - flag indicating matrix structure
453: Notes:
454: Recall that MatSetValues() uses 0-based row and column numbers
455: in Fortran as well as in C.
456: */
457: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
458: {
459: Mat A = *AA; /* Jacobian matrix */
460: AppCtx *appctx = (AppCtx *) ctx; /* user-defined application context */
461: PetscInt mstart = 0;
462: PetscInt mend = appctx->m;
464: PetscInt i, idx[3];
465: PetscScalar v[3], stwo = -2./(appctx->h*appctx->h), sone = -.5*stwo;
467: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
468: Compute entries for the locally owned part of the matrix
469: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
470: /*
471: Set matrix rows corresponding to boundary data
472: */
474: mstart = 0;
475: v[0] = 1.0;
476: MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
477: mstart++;
479: mend--;
480: v[0] = 1.0;
481: MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
483: /*
484: Set matrix rows corresponding to interior data. We construct the
485: matrix one row at a time.
486: */
487: v[0] = sone; v[1] = stwo; v[2] = sone;
488: for ( i=mstart; i<mend; i++ ) {
489: idx[0] = i-1; idx[1] = i; idx[2] = i+1;
490: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
491: }
493: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
494: Complete the matrix assembly process and set some options
495: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
496: /*
497: Assemble matrix, using the 2-step process:
498: MatAssemblyBegin(), MatAssemblyEnd()
499: Computations can be done while messages are in transition
500: by placing code between these two statements.
501: */
502: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
503: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
505: /*
506: Set flag to indicate that the Jacobian matrix retains an identical
507: nonzero structure throughout all timestepping iterations (although the
508: values of the entries change). Thus, we can save some work in setting
509: up the preconditioner (e.g., no need to redo symbolic factorization for
510: ILU/ICC preconditioners).
511: - If the nonzero structure of the matrix is different during
512: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
513: must be used instead. If you are unsure whether the matrix
514: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
515: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
516: believes your assertion and does not check the structure
517: of the matrix. If you erroneously claim that the structure
518: is the same when it actually is not, the new preconditioner
519: will not function correctly. Thus, use this optimization
520: feature with caution!
521: */
522: *str = SAME_NONZERO_PATTERN;
524: /*
525: Set and option to indicate that we will never add a new nonzero location
526: to the matrix. If we do, it will generate an error.
527: */
528: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
530: return 0;
531: }
532: /* --------------------------------------------------------------------- */
535: /*
536: Input Parameters:
537: ts - the TS context
538: t - current time
539: f - function
540: ctx - optional user-defined context, as set by TSetBCFunction()
541: */
542: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
543: {
544: AppCtx *appctx = (AppCtx *) ctx; /* user-defined application context */
546: PetscInt m = appctx->m;
547: PetscScalar *fa;
549: VecGetArray(f,&fa);
550: fa[0] = 0.0;
551: fa[m-1] = 1.0;
552: VecRestoreArray(f,&fa);
553: printf("t=%g\n",t);
554:
555: return 0;
556: }