Actual source code: ex5.c
petsc-3.3-p7 2013-05-11
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) = 1, u(t,1) = 1,
25: and the initial condition
26: u(0,x) = cos(6*pi*x) + 3*cos(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) * cos(6*pi*x) +
38: 3*exp(-4*pi*pi*t) * cos(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 just 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: */
78: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
79: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
80: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
81: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
85: int main(int argc,char **argv)
86: {
87: AppCtx appctx; /* user-defined application context */
88: TS ts; /* timestepping context */
89: Mat A; /* matrix data structure */
90: Vec u; /* approximate solution vector */
91: PetscReal time_total_max = 100.0; /* default max total time */
92: PetscInt time_steps_max = 100; /* default max timesteps */
93: PetscDraw draw; /* drawing context */
95: PetscInt steps,m;
96: PetscMPIInt size;
97: PetscBool flg;
98: PetscReal dt,ftime;
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Initialize program and set problem parameters
102: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:
104: PetscInitialize(&argc,&argv,(char*)0,help);
105: MPI_Comm_size(PETSC_COMM_WORLD,&size);
106: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
108: m = 60;
109: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
110: PetscOptionsHasName(PETSC_NULL,"-debug",&appctx.debug);
111: appctx.m = m;
112: appctx.h = 1.0/(m-1.0);
113: appctx.norm_2 = 0.0;
114: appctx.norm_max = 0.0;
115: PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");
117: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: Create vector data structures
119: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121: /*
122: Create vector data structures for approximate and exact solutions
123: */
124: VecCreateSeq(PETSC_COMM_SELF,m,&u);
125: VecDuplicate(u,&appctx.solution);
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Set up displays to show graphs of the solution and error
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
132: PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
133: PetscDrawSetDoubleBuffer(draw);
134: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
135: PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
136: PetscDrawSetDoubleBuffer(draw);
138: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: Create timestepping solver context
140: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: TSCreate(PETSC_COMM_SELF,&ts);
143: TSSetProblemType(ts,TS_LINEAR);
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Set optional user-defined monitoring routine
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: TSMonitorSet(ts,Monitor,&appctx,PETSC_NULL);
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Create matrix data structure; set matrix evaluation routine.
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: MatCreate(PETSC_COMM_SELF,&A);
157: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
158: MatSetFromOptions(A);
160: PetscOptionsHasName(PETSC_NULL,"-time_dependent_rhs",&flg);
161: if (flg) {
162: /*
163: For linear problems with a time-dependent f(u,t) in the equation
164: u_t = f(u,t), the user provides the discretized right-hand-side
165: as a time-dependent matrix.
166: */
167: TSSetRHSFunction(ts,PETSC_NULL,TSComputeRHSFunctionLinear,&appctx);
168: TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
169: } else {
170: /*
171: For linear problems with a time-independent f(u) in the equation
172: u_t = f(u), the user provides the discretized right-hand-side
173: as a matrix only once, and then sets a null matrix evaluation
174: routine.
175: */
176: MatStructure A_structure;
177: RHSMatrixHeat(ts,0.0,u,&A,&A,&A_structure,&appctx);
178: TSSetRHSFunction(ts,PETSC_NULL,TSComputeRHSFunctionLinear,&appctx);
179: TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
180: }
182: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183: Set solution vector and initial timestep
184: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186: dt = appctx.h*appctx.h/2.0;
187: TSSetInitialTimeStep(ts,0.0,dt);
188: TSSetSolution(ts,u);
190: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191: Customize timestepping solver:
192: - Set the solution method to be the Backward Euler method.
193: - Set timestepping duration info
194: Then set runtime options, which can override these defaults.
195: For example,
196: -ts_max_steps <maxsteps> -ts_final_time <maxtime>
197: to override the defaults set by TSSetDuration().
198: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200: TSSetDuration(ts,time_steps_max,time_total_max);
201: TSSetFromOptions(ts);
203: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204: Solve the problem
205: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207: /*
208: Evaluate initial conditions
209: */
210: InitialConditions(u,&appctx);
212: /*
213: Run the timestepping solver
214: */
215: TSSolve(ts,u,&ftime);
216: TSGetTimeStepNumber(ts,&steps);
218: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219: View timestepping solver info
220: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
222: PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %G, avg. error (max norm) = %G\n",
223: appctx.norm_2/steps,appctx.norm_max/steps);
224: TSView(ts,PETSC_VIEWER_STDOUT_SELF);
226: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227: Free work space. All PETSc objects should be destroyed when they
228: are no longer needed.
229: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231: TSDestroy(&ts);
232: MatDestroy(&A);
233: VecDestroy(&u);
234: PetscViewerDestroy(&appctx.viewer1);
235: PetscViewerDestroy(&appctx.viewer2);
236: VecDestroy(&appctx.solution);
238: /*
239: Always call PetscFinalize() before exiting a program. This routine
240: - finalizes the PETSc libraries as well as MPI
241: - provides summary and diagnostic information if certain runtime
242: options are chosen (e.g., -log_summary).
243: */
244: PetscFinalize();
245: return 0;
246: }
247: /* --------------------------------------------------------------------- */
250: /*
251: InitialConditions - Computes the solution at the initial time.
253: Input Parameter:
254: u - uninitialized solution vector (global)
255: appctx - user-defined application context
257: Output Parameter:
258: u - vector with solution at initial time (global)
259: */
260: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
261: {
262: PetscScalar *u_localptr,h = appctx->h;
263: PetscInt i;
266: /*
267: Get a pointer to vector data.
268: - For default PETSc vectors, VecGetArray() returns a pointer to
269: the data array. Otherwise, the routine is implementation dependent.
270: - You MUST call VecRestoreArray() when you no longer need access to
271: the array.
272: - Note that the Fortran interface to VecGetArray() differs from the
273: C version. See the users manual for details.
274: */
275: VecGetArray(u,&u_localptr);
277: /*
278: We initialize the solution array by simply writing the solution
279: directly into the array locations. Alternatively, we could use
280: VecSetValues() or VecSetValuesLocal().
281: */
282: for (i=0; i<appctx->m; i++) {
283: u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h);
284: }
286: /*
287: Restore vector
288: */
289: VecRestoreArray(u,&u_localptr);
291: /*
292: Print debugging information if desired
293: */
294: if (appctx->debug) {
295: printf("initial guess vector\n");
296: VecView(u,PETSC_VIEWER_STDOUT_SELF);
297: }
299: return 0;
300: }
301: /* --------------------------------------------------------------------- */
304: /*
305: ExactSolution - Computes the exact solution at a given time.
307: Input Parameters:
308: t - current time
309: solution - vector in which exact solution will be computed
310: appctx - user-defined application context
312: Output Parameter:
313: solution - vector with the newly computed exact solution
314: */
315: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
316: {
317: PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
318: PetscInt i;
321: /*
322: Get a pointer to vector data.
323: */
324: VecGetArray(solution,&s_localptr);
326: /*
327: Simply write the solution directly into the array locations.
328: Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
329: */
330: ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
331: sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h;
332: for (i=0; i<appctx->m; i++) {
333: s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2;
334: }
336: /*
337: Restore vector
338: */
339: VecRestoreArray(solution,&s_localptr);
340: return 0;
341: }
342: /* --------------------------------------------------------------------- */
345: /*
346: Monitor - User-provided routine to monitor the solution computed at
347: each timestep. This example plots the solution and computes the
348: error in two different norms.
350: Input Parameters:
351: ts - the timestep context
352: step - the count of the current step (with 0 meaning the
353: initial condition)
354: time - the current time
355: u - the solution at this timestep
356: ctx - the user-provided context for this monitoring routine.
357: In this case we use the application context which contains
358: information about the problem size, workspace and the exact
359: solution.
360: */
361: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
362: {
363: AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */
365: PetscReal norm_2,norm_max;
367: /*
368: View a graph of the current iterate
369: */
370: VecView(u,appctx->viewer2);
372: /*
373: Compute the exact solution
374: */
375: ExactSolution(time,appctx->solution,appctx);
377: /*
378: Print debugging information if desired
379: */
380: if (appctx->debug) {
381: printf("Computed solution vector\n");
382: VecView(u,PETSC_VIEWER_STDOUT_SELF);
383: printf("Exact solution vector\n");
384: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
385: }
387: /*
388: Compute the 2-norm and max-norm of the error
389: */
390: VecAXPY(appctx->solution,-1.0,u);
391: VecNorm(appctx->solution,NORM_2,&norm_2);
392: norm_2 = PetscSqrtReal(appctx->h)*norm_2;
393: VecNorm(appctx->solution,NORM_MAX,&norm_max);
395: PetscPrintf(PETSC_COMM_WORLD,"Timestep %D: time = %G, 2-norm error = %G, max norm error = %G\n",
396: step,time,norm_2,norm_max);
397: appctx->norm_2 += norm_2;
398: appctx->norm_max += norm_max;
400: /*
401: View a graph of the error
402: */
403: VecView(appctx->solution,appctx->viewer1);
405: /*
406: Print debugging information if desired
407: */
408: if (appctx->debug) {
409: printf("Error vector\n");
410: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
411: }
413: return 0;
414: }
415: /* --------------------------------------------------------------------- */
418: /*
419: RHSMatrixHeat - User-provided routine to compute the right-hand-side
420: matrix for the heat equation.
422: Input Parameters:
423: ts - the TS context
424: t - current time
425: global_in - global input vector
426: dummy - optional user-defined context, as set by TSetRHSJacobian()
428: Output Parameters:
429: AA - Jacobian matrix
430: BB - optionally different preconditioning matrix
431: str - flag indicating matrix structure
433: Notes:
434: Recall that MatSetValues() uses 0-based row and column numbers
435: in Fortran as well as in C.
436: */
437: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
438: {
439: Mat A = *AA; /* Jacobian matrix */
440: AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */
441: PetscInt mstart = 0;
442: PetscInt mend = appctx->m;
444: PetscInt i,idx[3];
445: PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;
447: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
448: Compute entries for the locally owned part of the matrix
449: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
450: /*
451: Set matrix rows corresponding to boundary data
452: */
454: mstart = 0;
455: v[0] = 1.0;
456: MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
457: mstart++;
459: mend--;
460: v[0] = 1.0;
461: MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
463: /*
464: Set matrix rows corresponding to interior data. We construct the
465: matrix one row at a time.
466: */
467: v[0] = sone; v[1] = stwo; v[2] = sone;
468: for (i=mstart; i<mend; i++) {
469: idx[0] = i-1; idx[1] = i; idx[2] = i+1;
470: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
471: }
473: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
474: Complete the matrix assembly process and set some options
475: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
476: /*
477: Assemble matrix, using the 2-step process:
478: MatAssemblyBegin(), MatAssemblyEnd()
479: Computations can be done while messages are in transition
480: by placing code between these two statements.
481: */
482: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
483: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
485: /*
486: Set flag to indicate that the Jacobian matrix retains an identical
487: nonzero structure throughout all timestepping iterations (although the
488: values of the entries change). Thus, we can save some work in setting
489: up the preconditioner (e.g., no need to redo symbolic factorization for
490: ILU/ICC preconditioners).
491: - If the nonzero structure of the matrix is different during
492: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
493: must be used instead. If you are unsure whether the matrix
494: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
495: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
496: believes your assertion and does not check the structure
497: of the matrix. If you erroneously claim that the structure
498: is the same when it actually is not, the new preconditioner
499: will not function correctly. Thus, use this optimization
500: feature with caution!
501: */
502: *str = SAME_NONZERO_PATTERN;
504: /*
505: Set and option to indicate that we will never add a new nonzero location
506: to the matrix. If we do, it will generate an error.
507: */
508: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
510: return 0;
511: }