Actual source code: ex3.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) = 0, u(t,1) = 0,
23: and the initial condition
24: u(0,x) = sin(6*pi*x) + 3*sin(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) * sin(6*pi*x) +
36: 3*exp(-4*pi*pi*t) * sin(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 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: */
58: #include <petscts.h>
59: #include <petscdraw.h>
61: /*
62: User-defined application context - contains data needed by the
63: application-provided call-back routines.
64: */
65: typedef struct {
66: Vec solution; /* global exact solution vector */
67: PetscInt m; /* total number of grid points */
68: PetscReal h; /* mesh width h = 1/(m-1) */
69: PetscBool debug; /* flag (1 indicates activation of debugging printouts) */
70: PetscViewer viewer1,viewer2; /* viewers for the solution and error */
71: PetscReal norm_2,norm_max; /* error norms */
72: } AppCtx;
74: /*
75: User-defined routines
76: */
77: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
78: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
79: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
80: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
84: int main(int argc,char **argv)
85: {
86: AppCtx appctx; /* user-defined application context */
87: TS ts; /* timestepping context */
88: Mat A; /* matrix data structure */
89: Vec u; /* approximate solution vector */
90: PetscReal time_total_max = 100.0; /* default max total time */
91: PetscInt time_steps_max = 100; /* default max timesteps */
92: PetscDraw draw; /* drawing context */
94: PetscInt steps,m;
95: PetscMPIInt size;
96: PetscReal dt;
97: PetscBool flg;
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Initialize program and set problem parameters
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: PetscInitialize(&argc,&argv,(char*)0,help);
104: MPI_Comm_size(PETSC_COMM_WORLD,&size);
105: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
107: m = 60;
108: PetscOptionsGetInt(NULL,"-m",&m,NULL);
109: PetscOptionsHasName(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;
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,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);
160: MatSetUp(A);
162: flg = PETSC_FALSE;
163: PetscOptionsGetBool(NULL,"-time_dependent_rhs",&flg,NULL);
164: if (flg) {
165: /*
166: For linear problems with a time-dependent f(u,t) in the equation
167: u_t = f(u,t), the user provides the discretized right-hand-side
168: as a time-dependent matrix.
169: */
170: TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
171: TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
172: } else {
173: /*
174: For linear problems with a time-independent f(u) in the equation
175: u_t = f(u), the user provides the discretized right-hand-side
176: as a matrix only once, and then sets the special Jacobian evaluation
177: routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian.
178: */
179: MatStructure A_structure;
180: RHSMatrixHeat(ts,0.0,u,&A,&A,&A_structure,&appctx);
181: TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
182: TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
183: }
185: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186: Set solution vector and initial timestep
187: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189: dt = appctx.h*appctx.h/2.0;
190: TSSetInitialTimeStep(ts,0.0,dt);
192: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193: Customize timestepping solver:
194: - Set the solution method to be the Backward Euler method.
195: - Set timestepping duration info
196: Then set runtime options, which can override these defaults.
197: For example,
198: -ts_max_steps <maxsteps> -ts_final_time <maxtime>
199: to override the defaults set by TSSetDuration().
200: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202: TSSetDuration(ts,time_steps_max,time_total_max);
203: TSSetFromOptions(ts);
205: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: Solve the problem
207: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209: /*
210: Evaluate initial conditions
211: */
212: InitialConditions(u,&appctx);
214: /*
215: Run the timestepping solver
216: */
217: TSSolve(ts,u);
218: TSGetTimeStepNumber(ts,&steps);
220: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221: View timestepping solver info
222: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224: PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %G, avg. error (max norm) = %G\n",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++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(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: PetscPrintf(PETSC_COMM_WORLD,"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;
318: 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);
330: 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++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
334: /*
335: Restore vector
336: */
337: VecRestoreArray(solution,&s_localptr);
338: return 0;
339: }
340: /* --------------------------------------------------------------------- */
343: /*
344: Monitor - User-provided routine to monitor the solution computed at
345: each timestep. This example plots the solution and computes the
346: error in two different norms.
348: This example also demonstrates changing the timestep via TSSetTimeStep().
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,dt,dttol;
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: PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n");
382: VecView(u,PETSC_VIEWER_STDOUT_SELF);
383: PetscPrintf(PETSC_COMM_SELF,"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: TSGetTimeStep(ts,&dt);
396: PetscPrintf(PETSC_COMM_WORLD,"Timestep %3D: step size = %-11g, time = %-11g, 2-norm error = %-11g, max norm error = %-11g\n",step,dt,time,norm_2,norm_max);
398: appctx->norm_2 += norm_2;
399: appctx->norm_max += norm_max;
401: dttol = .0001;
402: PetscOptionsGetReal(NULL,"-dttol",&dttol,NULL);
403: if (dt < dttol) {
404: dt *= .999;
405: TSSetTimeStep(ts,dt);
406: }
408: /*
409: View a graph of the error
410: */
411: VecView(appctx->solution,appctx->viewer1);
413: /*
414: Print debugging information if desired
415: */
416: if (appctx->debug) {
417: PetscPrintf(PETSC_COMM_SELF,"Error vector\n");
418: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
419: }
421: return 0;
422: }
423: /* --------------------------------------------------------------------- */
426: /*
427: RHSMatrixHeat - User-provided routine to compute the right-hand-side
428: matrix for the heat equation.
430: Input Parameters:
431: ts - the TS context
432: t - current time
433: global_in - global input vector
434: dummy - optional user-defined context, as set by TSetRHSJacobian()
436: Output Parameters:
437: AA - Jacobian matrix
438: BB - optionally different preconditioning matrix
439: str - flag indicating matrix structure
441: Notes:
442: Recall that MatSetValues() uses 0-based row and column numbers
443: in Fortran as well as in C.
444: */
445: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
446: {
447: Mat A = *AA; /* Jacobian matrix */
448: AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */
449: PetscInt mstart = 0;
450: PetscInt mend = appctx->m;
452: PetscInt i,idx[3];
453: PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;
455: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
456: Compute entries for the locally owned part of the matrix
457: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
458: /*
459: Set matrix rows corresponding to boundary data
460: */
462: mstart = 0;
463: v[0] = 1.0;
464: MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
465: mstart++;
467: mend--;
468: v[0] = 1.0;
469: MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
471: /*
472: Set matrix rows corresponding to interior data. We construct the
473: matrix one row at a time.
474: */
475: v[0] = sone; v[1] = stwo; v[2] = sone;
476: for (i=mstart; i<mend; i++) {
477: idx[0] = i-1; idx[1] = i; idx[2] = i+1;
478: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
479: }
481: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
482: Complete the matrix assembly process and set some options
483: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
484: /*
485: Assemble matrix, using the 2-step process:
486: MatAssemblyBegin(), MatAssemblyEnd()
487: Computations can be done while messages are in transition
488: by placing code between these two statements.
489: */
490: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
491: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
493: /*
494: Set flag to indicate that the Jacobian matrix retains an identical
495: nonzero structure throughout all timestepping iterations (although the
496: values of the entries change). Thus, we can save some work in setting
497: up the preconditioner (e.g., no need to redo symbolic factorization for
498: ILU/ICC preconditioners).
499: - If the nonzero structure of the matrix is different during
500: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
501: must be used instead. If you are unsure whether the matrix
502: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
503: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
504: believes your assertion and does not check the structure
505: of the matrix. If you erroneously claim that the structure
506: is the same when it actually is not, the new preconditioner
507: will not function correctly. Thus, use this optimization
508: feature with caution!
509: */
510: *str = SAME_NONZERO_PATTERN;
512: /*
513: Set and option to indicate that we will never add a new nonzero location
514: to the matrix. If we do, it will generate an error.
515: */
516: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
518: return 0;
519: }