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