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: petsc.h - base PETSc routines petscvec.h - vectors
54: petscsys.h - system routines 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: PetscTruth 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,ftime;
99: PetscTruth 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(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: PetscOptionsHasName(PETSC_NULL,"-time_dependent_rhs",&flg);
162: if (flg) {
163: /*
164: For linear problems with a time-dependent f(u,t) in the equation
165: u_t = f(u,t), the user provides the discretized right-hand-side
166: as a time-dependent matrix.
167: */
168: TSSetRHSMatrix(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,&A,&A,&A_structure,&appctx);
178: TSSetRHSMatrix(ts,A,A,PETSC_NULL,&appctx);
179: }
181: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182: Set solution vector and initial timestep
183: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185: dt = appctx.h*appctx.h/2.0;
186: TSSetInitialTimeStep(ts,0.0,dt);
187: TSSetSolution(ts,u);
189: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190: Customize timestepping solver:
191: - Set the solution method to be the Backward Euler method.
192: - Set timestepping duration info
193: Then set runtime options, which can override these defaults.
194: For example,
195: -ts_max_steps <maxsteps> -ts_max_time <maxtime>
196: to override the defaults set by TSSetDuration().
197: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
199: TSSetDuration(ts,time_steps_max,time_total_max);
200: TSSetFromOptions(ts);
202: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: Solve the problem
204: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206: /*
207: Evaluate initial conditions
208: */
209: InitialConditions(u,&appctx);
211: /*
212: Run the timestepping solver
213: */
214: TSStep(ts,&steps,&ftime);
216: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217: View timestepping solver info
218: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220: PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %G, avg. error (max norm) = %G\n",
221: appctx.norm_2/steps,appctx.norm_max/steps);
222: TSView(ts,PETSC_VIEWER_STDOUT_SELF);
224: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: Free work space. All PETSc objects should be destroyed when they
226: are no longer needed.
227: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229: TSDestroy(ts);
230: MatDestroy(A);
231: VecDestroy(u);
232: PetscViewerDestroy(appctx.viewer1);
233: PetscViewerDestroy(appctx.viewer2);
234: VecDestroy(appctx.solution);
236: /*
237: Always call PetscFinalize() before exiting a program. This routine
238: - finalizes the PETSc libraries as well as MPI
239: - provides summary and diagnostic information if certain runtime
240: options are chosen (e.g., -log_summary).
241: */
242: PetscFinalize();
243: return 0;
244: }
245: /* --------------------------------------------------------------------- */
248: /*
249: InitialConditions - Computes the solution at the initial time.
251: Input Parameter:
252: u - uninitialized solution vector (global)
253: appctx - user-defined application context
255: Output Parameter:
256: u - vector with solution at initial time (global)
257: */
258: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
259: {
260: PetscScalar *u_localptr,h = appctx->h;
262: PetscInt i;
264: /*
265: Get a pointer to vector data.
266: - For default PETSc vectors, VecGetArray() returns a pointer to
267: the data array. Otherwise, the routine is implementation dependent.
268: - You MUST call VecRestoreArray() when you no longer need access to
269: the array.
270: - Note that the Fortran interface to VecGetArray() differs from the
271: C version. See the users manual for details.
272: */
273: VecGetArray(u,&u_localptr);
275: /*
276: We initialize the solution array by simply writing the solution
277: directly into the array locations. Alternatively, we could use
278: VecSetValues() or VecSetValuesLocal().
279: */
280: for (i=0; i<appctx->m; i++) {
281: u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
282: }
284: /*
285: Restore vector
286: */
287: VecRestoreArray(u,&u_localptr);
289: /*
290: Print debugging information if desired
291: */
292: if (appctx->debug) {
293: printf("initial guess vector\n");
294: VecView(u,PETSC_VIEWER_STDOUT_SELF);
295: }
297: return 0;
298: }
299: /* --------------------------------------------------------------------- */
302: /*
303: ExactSolution - Computes the exact solution at a given time.
305: Input Parameters:
306: t - current time
307: solution - vector in which exact solution will be computed
308: appctx - user-defined application context
310: Output Parameter:
311: solution - vector with the newly computed exact solution
312: */
313: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
314: {
315: PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
317: PetscInt i;
319: /*
320: Get a pointer to vector data.
321: */
322: VecGetArray(solution,&s_localptr);
324: /*
325: Simply write the solution directly into the array locations.
326: Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
327: */
328: ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc);
329: 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++) {
332: s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
333: }
335: /*
336: Restore vector
337: */
338: VecRestoreArray(solution,&s_localptr);
339: return 0;
340: }
341: /* --------------------------------------------------------------------- */
344: /*
345: Monitor - User-provided routine to monitor the solution computed at
346: each timestep. This example plots the solution and computes the
347: error in two different norms.
349: This example also demonstrates changing the timestep via TSSetTimeStep().
351: Input Parameters:
352: ts - the timestep context
353: step - the count of the current step (with 0 meaning the
354: initial condition)
355: time - the current time
356: u - the solution at this timestep
357: ctx - the user-provided context for this monitoring routine.
358: In this case we use the application context which contains
359: information about the problem size, workspace and the exact
360: solution.
361: */
362: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
363: {
364: AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */
366: 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: 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 = sqrt(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",
397: step,dt,time,norm_2,norm_max);
398: appctx->norm_2 += norm_2;
399: appctx->norm_max += norm_max;
401: dttol = .0001;
402: PetscOptionsGetReal(PETSC_NULL,"-dttol",&dttol,PETSC_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: printf("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,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);
518: return 0;
519: }
520: /* --------------------------------------------------------------------- */
523: /*
524: Input Parameters:
525: ts - the TS context
526: t - current time
527: f - function
528: ctx - optional user-defined context, as set by TSetBCFunction()
529: */
530: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
531: {
532: AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */
533: PetscErrorCode ierr,m = appctx->m;
534: PetscScalar *fa;
536: VecGetArray(f,&fa);
537: fa[0] = 0.0;
538: fa[m-1] = 0.0;
539: VecRestoreArray(f,&fa);
540: printf("t=%g\n",t);
541:
542: return 0;
543: }