Actual source code: ex5.c
petsc-3.4.5 2014-06-29
2: static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\
3: Input parameters include:\n\
4: -m <points>, where <points> = number of grid points\n\
5: -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
6: -debug : Activate debugging printouts\n\
7: -nox : Deactivate x-window graphics\n\n";
9: /*
10: Concepts: TS^time-dependent linear problems
11: Concepts: TS^heat equation
12: Concepts: TS^diffusion equation
13: Processors: 1
14: */
16: /* ------------------------------------------------------------------------
18: This program solves the one-dimensional heat equation (also called the
19: diffusion equation),
20: u_t = u_xx,
21: on the domain 0 <= x <= 1, with the boundary conditions
22: u(t,0) = 1, u(t,1) = 1,
23: and the initial condition
24: u(0,x) = cos(6*pi*x) + 3*cos(2*pi*x).
25: This is a linear, second-order, parabolic equation.
27: We discretize the right-hand side using finite differences with
28: uniform grid spacing h:
29: u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
30: We then demonstrate time evolution using the various TS methods by
31: running the program via
32: ex3 -ts_type <timestepping solver>
34: We compare the approximate solution with the exact solution, given by
35: u_exact(x,t) = exp(-36*pi*pi*t) * cos(6*pi*x) +
36: 3*exp(-4*pi*pi*t) * cos(2*pi*x)
38: Notes:
39: This code demonstrates the TS solver interface to two variants of
40: linear problems, u_t = f(u,t), namely
41: - time-dependent f: f(u,t) is a function of t
42: - time-independent f: f(u,t) is simply just f(u)
44: The parallel version of this code is ts/examples/tutorials/ex4.c
46: ------------------------------------------------------------------------- */
48: /*
49: Include "petscts.h" so that we can use TS solvers. Note that this file
50: automatically includes:
51: petscsys.h - base PETSc routines petscvec.h - vectors
52: petscmat.h - matrices
53: petscis.h - index sets petscksp.h - Krylov subspace methods
54: petscviewer.h - viewers petscpc.h - preconditioners
55: petscksp.h - linear solvers petscsnes.h - nonlinear solvers
56: */
57: #include <petscts.h>
58: #include <petscdraw.h>
60: /*
61: User-defined application context - contains data needed by the
62: application-provided call-back routines.
63: */
64: typedef struct {
65: Vec solution; /* global exact solution vector */
66: PetscInt m; /* total number of grid points */
67: PetscReal h; /* mesh width h = 1/(m-1) */
68: PetscBool debug; /* flag (1 indicates activation of debugging printouts) */
69: PetscViewer viewer1,viewer2; /* viewers for the solution and error */
70: PetscReal norm_2,norm_max; /* error norms */
71: } AppCtx;
73: /*
74: User-defined routines
75: */
76: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
77: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
78: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
79: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
83: int main(int argc,char **argv)
84: {
85: AppCtx appctx; /* user-defined application context */
86: TS ts; /* timestepping context */
87: Mat A; /* matrix data structure */
88: Vec u; /* approximate solution vector */
89: PetscReal time_total_max = 100.0; /* default max total time */
90: PetscInt time_steps_max = 100; /* default max timesteps */
91: PetscDraw draw; /* drawing context */
93: PetscInt steps,m;
94: PetscMPIInt size;
95: PetscBool flg;
96: PetscReal dt,ftime;
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Initialize program and set problem parameters
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: PetscInitialize(&argc,&argv,(char*)0,help);
103: MPI_Comm_size(PETSC_COMM_WORLD,&size);
104: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
106: m = 60;
107: PetscOptionsGetInt(NULL,"-m",&m,NULL);
108: PetscOptionsHasName(NULL,"-debug",&appctx.debug);
109: appctx.m = m;
110: appctx.h = 1.0/(m-1.0);
111: appctx.norm_2 = 0.0;
112: appctx.norm_max = 0.0;
114: PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Create vector data structures
118: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120: /*
121: Create vector data structures for approximate and exact solutions
122: */
123: VecCreateSeq(PETSC_COMM_SELF,m,&u);
124: VecDuplicate(u,&appctx.solution);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Set up displays to show graphs of the solution and error
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
131: PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
132: PetscDrawSetDoubleBuffer(draw);
133: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
134: PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
135: PetscDrawSetDoubleBuffer(draw);
137: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: Create timestepping solver context
139: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141: TSCreate(PETSC_COMM_SELF,&ts);
142: TSSetProblemType(ts,TS_LINEAR);
144: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145: Set optional user-defined monitoring routine
146: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148: TSMonitorSet(ts,Monitor,&appctx,NULL);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Create matrix data structure; set matrix evaluation routine.
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: MatCreate(PETSC_COMM_SELF,&A);
156: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
157: MatSetFromOptions(A);
158: MatSetUp(A);
160: PetscOptionsHasName(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,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,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);
216: TSGetSolveTime(ts,&ftime);
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;
264: 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++) u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h);
285: /*
286: Restore vector
287: */
288: VecRestoreArray(u,&u_localptr);
290: /*
291: Print debugging information if desired
292: */
293: if (appctx->debug) {
294: printf("initial guess vector\n");
295: VecView(u,PETSC_VIEWER_STDOUT_SELF);
296: }
298: return 0;
299: }
300: /* --------------------------------------------------------------------- */
303: /*
304: ExactSolution - Computes the exact solution at a given time.
306: Input Parameters:
307: t - current time
308: solution - vector in which exact solution will be computed
309: appctx - user-defined application context
311: Output Parameter:
312: solution - vector with the newly computed exact solution
313: */
314: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
315: {
316: PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
317: PetscInt i;
320: /*
321: Get a pointer to vector data.
322: */
323: VecGetArray(solution,&s_localptr);
325: /*
326: Simply write the solution directly into the array locations.
327: Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
328: */
329: ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
330: sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h;
331: for (i=0; i<appctx->m; i++) s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2;
333: /*
334: Restore vector
335: */
336: VecRestoreArray(solution,&s_localptr);
337: return 0;
338: }
339: /* --------------------------------------------------------------------- */
342: /*
343: Monitor - User-provided routine to monitor the solution computed at
344: each timestep. This example plots the solution and computes the
345: error in two different norms.
347: Input Parameters:
348: ts - the timestep context
349: step - the count of the current step (with 0 meaning the
350: initial condition)
351: time - the current time
352: u - the solution at this timestep
353: ctx - the user-provided context for this monitoring routine.
354: In this case we use the application context which contains
355: information about the problem size, workspace and the exact
356: solution.
357: */
358: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
359: {
360: AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */
362: PetscReal norm_2,norm_max;
364: /*
365: View a graph of the current iterate
366: */
367: VecView(u,appctx->viewer2);
369: /*
370: Compute the exact solution
371: */
372: ExactSolution(time,appctx->solution,appctx);
374: /*
375: Print debugging information if desired
376: */
377: if (appctx->debug) {
378: printf("Computed solution vector\n");
379: VecView(u,PETSC_VIEWER_STDOUT_SELF);
380: printf("Exact solution vector\n");
381: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
382: }
384: /*
385: Compute the 2-norm and max-norm of the error
386: */
387: VecAXPY(appctx->solution,-1.0,u);
388: VecNorm(appctx->solution,NORM_2,&norm_2);
389: norm_2 = PetscSqrtReal(appctx->h)*norm_2;
390: VecNorm(appctx->solution,NORM_MAX,&norm_max);
392: PetscPrintf(PETSC_COMM_WORLD,"Timestep %D: time = %G, 2-norm error = %G, max norm error = %G\n",
393: step,time,norm_2,norm_max);
394: appctx->norm_2 += norm_2;
395: appctx->norm_max += norm_max;
397: /*
398: View a graph of the error
399: */
400: VecView(appctx->solution,appctx->viewer1);
402: /*
403: Print debugging information if desired
404: */
405: if (appctx->debug) {
406: printf("Error vector\n");
407: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
408: }
410: return 0;
411: }
412: /* --------------------------------------------------------------------- */
415: /*
416: RHSMatrixHeat - User-provided routine to compute the right-hand-side
417: matrix for the heat equation.
419: Input Parameters:
420: ts - the TS context
421: t - current time
422: global_in - global input vector
423: dummy - optional user-defined context, as set by TSetRHSJacobian()
425: Output Parameters:
426: AA - Jacobian matrix
427: BB - optionally different preconditioning matrix
428: str - flag indicating matrix structure
430: Notes:
431: Recall that MatSetValues() uses 0-based row and column numbers
432: in Fortran as well as in C.
433: */
434: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
435: {
436: Mat A = *AA; /* Jacobian matrix */
437: AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */
438: PetscInt mstart = 0;
439: PetscInt mend = appctx->m;
441: PetscInt i,idx[3];
442: PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;
444: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
445: Compute entries for the locally owned part of the matrix
446: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
447: /*
448: Set matrix rows corresponding to boundary data
449: */
451: mstart = 0;
452: v[0] = 1.0;
453: MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
454: mstart++;
456: mend--;
457: v[0] = 1.0;
458: MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
460: /*
461: Set matrix rows corresponding to interior data. We construct the
462: matrix one row at a time.
463: */
464: v[0] = sone; v[1] = stwo; v[2] = sone;
465: for (i=mstart; i<mend; i++) {
466: idx[0] = i-1; idx[1] = i; idx[2] = i+1;
467: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
468: }
470: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
471: Complete the matrix assembly process and set some options
472: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
473: /*
474: Assemble matrix, using the 2-step process:
475: MatAssemblyBegin(), MatAssemblyEnd()
476: Computations can be done while messages are in transition
477: by placing code between these two statements.
478: */
479: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
480: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
482: /*
483: Set flag to indicate that the Jacobian matrix retains an identical
484: nonzero structure throughout all timestepping iterations (although the
485: values of the entries change). Thus, we can save some work in setting
486: up the preconditioner (e.g., no need to redo symbolic factorization for
487: ILU/ICC preconditioners).
488: - If the nonzero structure of the matrix is different during
489: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
490: must be used instead. If you are unsure whether the matrix
491: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
492: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
493: believes your assertion and does not check the structure
494: of the matrix. If you erroneously claim that the structure
495: is the same when it actually is not, the new preconditioner
496: will not function correctly. Thus, use this optimization
497: feature with caution!
498: */
499: *str = SAME_NONZERO_PATTERN;
501: /*
502: Set and option to indicate that we will never add a new nonzero location
503: to the matrix. If we do, it will generate an error.
504: */
505: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
507: return 0;
508: }