Actual source code: ex3.c
petsc-3.10.5 2019-03-28
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: -use_ifunc : Use IFunction/IJacobian interface\n\
7: -debug : Activate debugging printouts\n\
8: -nox : Deactivate x-window graphics\n\n";
10: /*
11: Concepts: TS^time-dependent linear problems
12: Concepts: TS^heat equation
13: Concepts: TS^diffusion equation
14: Processors: 1
15: */
17: /* ------------------------------------------------------------------------
19: This program solves the one-dimensional heat equation (also called the
20: diffusion equation),
21: u_t = u_xx,
22: on the domain 0 <= x <= 1, with the boundary conditions
23: u(t,0) = 0, u(t,1) = 0,
24: and the initial condition
25: u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
26: This is a linear, second-order, parabolic equation.
28: We discretize the right-hand side using finite differences with
29: uniform grid spacing h:
30: u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
31: We then demonstrate time evolution using the various TS methods by
32: running the program via
33: ex3 -ts_type <timestepping solver>
35: We compare the approximate solution with the exact solution, given by
36: u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
37: 3*exp(-4*pi*pi*t) * sin(2*pi*x)
39: Notes:
40: This code demonstrates the TS solver interface to two variants of
41: linear problems, u_t = f(u,t), namely
42: - time-dependent f: f(u,t) is a function of t
43: - time-independent f: f(u,t) is simply f(u)
45: The parallel version of this code is ts/examples/tutorials/ex4.c
47: ------------------------------------------------------------------------- */
49: /*
50: Include "petscts.h" so that we can use TS solvers. Note that this file
51: automatically includes:
52: petscsys.h - base PETSc routines petscvec.h - vectors
53: petscmat.h - matrices
54: petscis.h - index sets petscksp.h - Krylov subspace methods
55: petscviewer.h - viewers petscpc.h - preconditioners
56: petscksp.h - linear solvers petscsnes.h - nonlinear solvers
57: */
59: #include <petscts.h>
60: #include <petscdraw.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: Mat A; /* RHS mat, used with IFunction interface */
74: PetscReal oshift; /* old shift applied, prevent to recompute the IJacobian */
75: } AppCtx;
77: /*
78: User-defined routines
79: */
80: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
81: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
82: extern PetscErrorCode IFunctionHeat(TS,PetscReal,Vec,Vec,Vec,void*);
83: extern PetscErrorCode IJacobianHeat(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
84: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
85: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
87: int main(int argc,char **argv)
88: {
89: AppCtx appctx; /* user-defined application context */
90: TS ts; /* timestepping context */
91: Mat A; /* matrix data structure */
92: Vec u; /* approximate solution vector */
93: PetscReal time_total_max = 100.0; /* default max total time */
94: PetscInt time_steps_max = 100; /* default max timesteps */
95: PetscDraw draw; /* drawing context */
97: PetscInt steps,m;
98: PetscMPIInt size;
99: PetscReal dt;
100: PetscBool flg;
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Initialize program and set problem parameters
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
107: MPI_Comm_size(PETSC_COMM_WORLD,&size);
108: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
110: m = 60;
111: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
112: PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);
114: appctx.m = m;
115: appctx.h = 1.0/(m-1.0);
116: appctx.norm_2 = 0.0;
117: appctx.norm_max = 0.0;
119: PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");
121: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: Create vector data structures
123: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: /*
126: Create vector data structures for approximate and exact solutions
127: */
128: VecCreateSeq(PETSC_COMM_SELF,m,&u);
129: VecDuplicate(u,&appctx.solution);
131: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132: Set up displays to show graphs of the solution and error
133: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
136: PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
137: PetscDrawSetDoubleBuffer(draw);
138: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
139: PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
140: PetscDrawSetDoubleBuffer(draw);
142: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: Create timestepping solver context
144: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146: TSCreate(PETSC_COMM_SELF,&ts);
147: TSSetProblemType(ts,TS_LINEAR);
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Set optional user-defined monitoring routine
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: TSMonitorSet(ts,Monitor,&appctx,NULL);
155: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: Create matrix data structure; set matrix evaluation routine.
158: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160: MatCreate(PETSC_COMM_SELF,&A);
161: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
162: MatSetFromOptions(A);
163: MatSetUp(A);
165: flg = PETSC_FALSE;
166: PetscOptionsGetBool(NULL,NULL,"-use_ifunc",&flg,NULL);
167: if (!flg) {
168: appctx.A = NULL;
169: PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL);
170: if (flg) {
171: /*
172: For linear problems with a time-dependent f(u,t) in the equation
173: u_t = f(u,t), the user provides the discretized right-hand-side
174: as a time-dependent matrix.
175: */
176: TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
177: TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
178: } else {
179: /*
180: For linear problems with a time-independent f(u) in the equation
181: u_t = f(u), the user provides the discretized right-hand-side
182: as a matrix only once, and then sets the special Jacobian evaluation
183: routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian.
184: */
185: RHSMatrixHeat(ts,0.0,u,A,A,&appctx);
186: TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
187: TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
188: }
189: } else {
190: Mat J;
192: RHSMatrixHeat(ts,0.0,u,A,A,&appctx);
193: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&J);
194: TSSetIFunction(ts,NULL,IFunctionHeat,&appctx);
195: TSSetIJacobian(ts,J,J,IJacobianHeat,&appctx);
196: MatDestroy(&J);
198: PetscObjectReference((PetscObject)A);
199: appctx.A = A;
200: appctx.oshift = PETSC_MIN_REAL;
201: }
202: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: Set solution vector and initial timestep
204: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206: dt = appctx.h*appctx.h/2.0;
207: TSSetTimeStep(ts,dt);
209: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210: Customize timestepping solver:
211: - Set the solution method to be the Backward Euler method.
212: - Set timestepping duration info
213: Then set runtime options, which can override these defaults.
214: For example,
215: -ts_max_steps <maxsteps> -ts_final_time <maxtime>
216: to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
217: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219: TSSetMaxSteps(ts,time_steps_max);
220: TSSetMaxTime(ts,time_total_max);
221: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
222: TSSetFromOptions(ts);
224: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: Solve the problem
226: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228: /*
229: Evaluate initial conditions
230: */
231: InitialConditions(u,&appctx);
233: /*
234: Run the timestepping solver
235: */
236: TSSolve(ts,u);
237: TSGetStepNumber(ts,&steps);
239: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
240: View timestepping solver info
241: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243: PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %g, avg. error (max norm) = %g\n",(double)(appctx.norm_2/steps),(double)(appctx.norm_max/steps));
244: TSView(ts,PETSC_VIEWER_STDOUT_SELF);
246: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247: Free work space. All PETSc objects should be destroyed when they
248: are no longer needed.
249: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
251: TSDestroy(&ts);
252: MatDestroy(&A);
253: VecDestroy(&u);
254: PetscViewerDestroy(&appctx.viewer1);
255: PetscViewerDestroy(&appctx.viewer2);
256: VecDestroy(&appctx.solution);
257: MatDestroy(&appctx.A);
259: /*
260: Always call PetscFinalize() before exiting a program. This routine
261: - finalizes the PETSc libraries as well as MPI
262: - provides summary and diagnostic information if certain runtime
263: options are chosen (e.g., -log_view).
264: */
265: PetscFinalize();
266: return ierr;
267: }
268: /* --------------------------------------------------------------------- */
269: /*
270: InitialConditions - Computes the solution at the initial time.
272: Input Parameter:
273: u - uninitialized solution vector (global)
274: appctx - user-defined application context
276: Output Parameter:
277: u - vector with solution at initial time (global)
278: */
279: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
280: {
281: PetscScalar *u_localptr,h = appctx->h;
283: PetscInt i;
285: /*
286: Get a pointer to vector data.
287: - For default PETSc vectors, VecGetArray() returns a pointer to
288: the data array. Otherwise, the routine is implementation dependent.
289: - You MUST call VecRestoreArray() when you no longer need access to
290: the array.
291: - Note that the Fortran interface to VecGetArray() differs from the
292: C version. See the users manual for details.
293: */
294: VecGetArray(u,&u_localptr);
296: /*
297: We initialize the solution array by simply writing the solution
298: directly into the array locations. Alternatively, we could use
299: VecSetValues() or VecSetValuesLocal().
300: */
301: for (i=0; i<appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
303: /*
304: Restore vector
305: */
306: VecRestoreArray(u,&u_localptr);
308: /*
309: Print debugging information if desired
310: */
311: if (appctx->debug) {
312: PetscPrintf(PETSC_COMM_WORLD,"Initial guess vector\n");
313: VecView(u,PETSC_VIEWER_STDOUT_SELF);
314: }
316: return 0;
317: }
318: /* --------------------------------------------------------------------- */
319: /*
320: ExactSolution - Computes the exact solution at a given time.
322: Input Parameters:
323: t - current time
324: solution - vector in which exact solution will be computed
325: appctx - user-defined application context
327: Output Parameter:
328: solution - vector with the newly computed exact solution
329: */
330: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
331: {
332: PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
334: PetscInt i;
336: /*
337: Get a pointer to vector data.
338: */
339: VecGetArray(solution,&s_localptr);
341: /*
342: Simply write the solution directly into the array locations.
343: Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
344: */
345: ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc);
346: ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
347: sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h;
348: for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
350: /*
351: Restore vector
352: */
353: VecRestoreArray(solution,&s_localptr);
354: return 0;
355: }
356: /* --------------------------------------------------------------------- */
357: /*
358: Monitor - User-provided routine to monitor the solution computed at
359: each timestep. This example plots the solution and computes the
360: error in two different norms.
362: This example also demonstrates changing the timestep via TSSetTimeStep().
364: Input Parameters:
365: ts - the timestep context
366: step - the count of the current step (with 0 meaning the
367: initial condition)
368: time - the current time
369: u - the solution at this timestep
370: ctx - the user-provided context for this monitoring routine.
371: In this case we use the application context which contains
372: information about the problem size, workspace and the exact
373: solution.
374: */
375: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
376: {
377: AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */
379: PetscReal norm_2,norm_max,dt,dttol;
381: /*
382: View a graph of the current iterate
383: */
384: VecView(u,appctx->viewer2);
386: /*
387: Compute the exact solution
388: */
389: ExactSolution(time,appctx->solution,appctx);
391: /*
392: Print debugging information if desired
393: */
394: if (appctx->debug) {
395: PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n");
396: VecView(u,PETSC_VIEWER_STDOUT_SELF);
397: PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");
398: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
399: }
401: /*
402: Compute the 2-norm and max-norm of the error
403: */
404: VecAXPY(appctx->solution,-1.0,u);
405: VecNorm(appctx->solution,NORM_2,&norm_2);
406: norm_2 = PetscSqrtReal(appctx->h)*norm_2;
407: VecNorm(appctx->solution,NORM_MAX,&norm_max);
408: if (norm_2 < 1e-14) norm_2 = 0;
409: if (norm_max < 1e-14) norm_max = 0;
411: TSGetTimeStep(ts,&dt);
412: PetscPrintf(PETSC_COMM_WORLD,"Timestep %3D: step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)dt,(double)time,(double)norm_2,(double)norm_max);
414: appctx->norm_2 += norm_2;
415: appctx->norm_max += norm_max;
417: dttol = .0001;
418: PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,NULL);
419: if (dt < dttol) {
420: dt *= .999;
421: TSSetTimeStep(ts,dt);
422: }
424: /*
425: View a graph of the error
426: */
427: VecView(appctx->solution,appctx->viewer1);
429: /*
430: Print debugging information if desired
431: */
432: if (appctx->debug) {
433: PetscPrintf(PETSC_COMM_SELF,"Error vector\n");
434: VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
435: }
437: return 0;
438: }
439: /* --------------------------------------------------------------------- */
440: /*
441: RHSMatrixHeat - User-provided routine to compute the right-hand-side
442: matrix for the heat equation.
444: Input Parameters:
445: ts - the TS context
446: t - current time
447: global_in - global input vector
448: dummy - optional user-defined context, as set by TSetRHSJacobian()
450: Output Parameters:
451: AA - Jacobian matrix
452: BB - optionally different preconditioning matrix
453: str - flag indicating matrix structure
455: Notes:
456: Recall that MatSetValues() uses 0-based row and column numbers
457: in Fortran as well as in C.
458: */
459: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
460: {
461: Mat A = AA; /* Jacobian matrix */
462: AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */
463: PetscInt mstart = 0;
464: PetscInt mend = appctx->m;
466: PetscInt i,idx[3];
467: PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;
469: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
470: Compute entries for the locally owned part of the matrix
471: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
472: /*
473: Set matrix rows corresponding to boundary data
474: */
476: mstart = 0;
477: v[0] = 1.0;
478: MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
479: mstart++;
481: mend--;
482: v[0] = 1.0;
483: MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
485: /*
486: Set matrix rows corresponding to interior data. We construct the
487: matrix one row at a time.
488: */
489: v[0] = sone; v[1] = stwo; v[2] = sone;
490: for (i=mstart; i<mend; i++) {
491: idx[0] = i-1; idx[1] = i; idx[2] = i+1;
492: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
493: }
495: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
496: Complete the matrix assembly process and set some options
497: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
498: /*
499: Assemble matrix, using the 2-step process:
500: MatAssemblyBegin(), MatAssemblyEnd()
501: Computations can be done while messages are in transition
502: by placing code between these two statements.
503: */
504: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
505: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
507: /*
508: Set and option to indicate that we will never add a new nonzero location
509: to the matrix. If we do, it will generate an error.
510: */
511: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
513: return 0;
514: }
516: PetscErrorCode IFunctionHeat(TS ts,PetscReal t,Vec X,Vec Xdot,Vec r,void *ctx)
517: {
518: AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */
521: MatMult(appctx->A,X,r);
522: VecAYPX(r,-1.0,Xdot);
523: return 0;
524: }
526: PetscErrorCode IJacobianHeat(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal s,Mat A,Mat B,void *ctx)
527: {
528: AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */
531: if (appctx->oshift == s) return 0;
532: MatCopy(appctx->A,A,SAME_NONZERO_PATTERN);
533: MatScale(A,-1);
534: MatShift(A,s);
535: MatCopy(A,B,SAME_NONZERO_PATTERN);
536: appctx->oshift = s;
537: return 0;
538: }
540: /*TEST
542: test:
543: args: -nox -ts_type ssp -ts_dt 0.0005
545: test:
546: suffix: 2
547: args: -nox -ts_type ssp -ts_dt 0.0005 -time_dependent_rhs 1
549: test:
550: suffix: 3
551: args: -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason
552: filter: sed "s/ATOL/RTOL/g"
553: requires: !single
555: test:
556: suffix: 4
557: args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason
558: filter: sed "s/ATOL/RTOL/g"
560: test:
561: suffix: 5
562: args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs
563: filter: sed "s/ATOL/RTOL/g"
565: test:
566: requires: !single
567: suffix: pod_guess
568: args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason
570: test:
571: requires: !single
572: suffix: pod_guess_Ainner
573: args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -ksp_guess_pod_Ainner -pc_type none -ksp_converged_reason
575: test:
576: requires: !single
577: suffix: fischer_guess
578: args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason
580: test:
581: requires: !single
582: suffix: fischer_guess_2
583: args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 2,10 -pc_type none -ksp_converged_reason
584: TEST*/