Actual source code: ex5.c


  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: /* ------------------------------------------------------------------------

 11:    This program solves the one-dimensional heat equation (also called the
 12:    diffusion equation),
 13:        u_t = u_xx,
 14:    on the domain 0 <= x <= 1, with the boundary conditions
 15:        u(t,0) = 1, u(t,1) = 1,
 16:    and the initial condition
 17:        u(0,x) = cos(6*pi*x) + 3*cos(2*pi*x).
 18:    This is a linear, second-order, parabolic equation.

 20:    We discretize the right-hand side using finite differences with
 21:    uniform grid spacing h:
 22:        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
 23:    We then demonstrate time evolution using the various TS methods by
 24:    running the program via
 25:        ex3 -ts_type <timestepping solver>

 27:    We compare the approximate solution with the exact solution, given by
 28:        u_exact(x,t) = exp(-36*pi*pi*t) * cos(6*pi*x) +
 29:                       3*exp(-4*pi*pi*t) * cos(2*pi*x)

 31:    Notes:
 32:    This code demonstrates the TS solver interface to two variants of
 33:    linear problems, u_t = f(u,t), namely
 34:      - time-dependent f:   f(u,t) is a function of t
 35:      - time-independent f: f(u,t) is simply just f(u)

 37:     The parallel version of this code is ts/tutorials/ex4.c

 39:   ------------------------------------------------------------------------- */

 41: /*
 42:    Include "petscts.h" so that we can use TS solvers.  Note that this file
 43:    automatically includes:
 44:      petscsys.h       - base PETSc routines   petscvec.h  - vectors
 45:      petscmat.h  - matrices
 46:      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
 47:      petscviewer.h - viewers               petscpc.h   - preconditioners
 48:      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
 49: */
 50: #include <petscts.h>
 51: #include <petscdraw.h>

 53: /*
 54:    User-defined application context - contains data needed by the
 55:    application-provided call-back routines.
 56: */
 57: typedef struct {
 58:   Vec         solution;          /* global exact solution vector */
 59:   PetscInt    m;                      /* total number of grid points */
 60:   PetscReal   h;                 /* mesh width h = 1/(m-1) */
 61:   PetscBool   debug;             /* flag (1 indicates activation of debugging printouts) */
 62:   PetscViewer viewer1,viewer2;  /* viewers for the solution and error */
 63:   PetscReal   norm_2,norm_max;  /* error norms */
 64: } AppCtx;

 66: /*
 67:    User-defined routines
 68: */
 69: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
 70: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
 71: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 72: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);

 74: int main(int argc,char **argv)
 75: {
 76:   AppCtx         appctx;                 /* user-defined application context */
 77:   TS             ts;                     /* timestepping context */
 78:   Mat            A;                      /* matrix data structure */
 79:   Vec            u;                      /* approximate solution vector */
 80:   PetscReal      time_total_max = 100.0; /* default max total time */
 81:   PetscInt       time_steps_max = 100;   /* default max timesteps */
 82:   PetscDraw      draw;                   /* drawing context */
 83:   PetscInt       steps,m;
 84:   PetscMPIInt    size;
 85:   PetscBool      flg;
 86:   PetscReal      dt,ftime;

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:      Initialize program and set problem parameters
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   PetscInitialize(&argc,&argv,(char*)0,help);
 93:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 96:   m               = 60;
 97:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 98:   PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);
 99:   appctx.m        = m;
100:   appctx.h        = 1.0/(m-1.0);
101:   appctx.norm_2   = 0.0;
102:   appctx.norm_max = 0.0;

104:   PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");

106:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107:      Create vector data structures
108:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

110:   /*
111:      Create vector data structures for approximate and exact solutions
112:   */
113:   VecCreateSeq(PETSC_COMM_SELF,m,&u);
114:   VecDuplicate(u,&appctx.solution);

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:      Set up displays to show graphs of the solution and error
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

120:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
121:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
122:   PetscDrawSetDoubleBuffer(draw);
123:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
124:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
125:   PetscDrawSetDoubleBuffer(draw);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:      Create timestepping solver context
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

131:   TSCreate(PETSC_COMM_SELF,&ts);
132:   TSSetProblemType(ts,TS_LINEAR);

134:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135:      Set optional user-defined monitoring routine
136:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

138:   TSMonitorSet(ts,Monitor,&appctx,NULL);

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

142:      Create matrix data structure; set matrix evaluation routine.
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

145:   MatCreate(PETSC_COMM_SELF,&A);
146:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
147:   MatSetFromOptions(A);
148:   MatSetUp(A);

150:   PetscOptionsHasName(NULL,NULL,"-time_dependent_rhs",&flg);
151:   if (flg) {
152:     /*
153:        For linear problems with a time-dependent f(u,t) in the equation
154:        u_t = f(u,t), the user provides the discretized right-hand-side
155:        as a time-dependent matrix.
156:     */
157:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
158:     TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
159:   } else {
160:     /*
161:        For linear problems with a time-independent f(u) in the equation
162:        u_t = f(u), the user provides the discretized right-hand-side
163:        as a matrix only once, and then sets a null matrix evaluation
164:        routine.
165:     */
166:     RHSMatrixHeat(ts,0.0,u,A,A,&appctx);
167:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
168:     TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
169:   }

171:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172:      Set solution vector and initial timestep
173:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

175:   dt   = appctx.h*appctx.h/2.0;
176:   TSSetTimeStep(ts,dt);
177:   TSSetSolution(ts,u);

179:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180:      Customize timestepping solver:
181:        - Set the solution method to be the Backward Euler method.
182:        - Set timestepping duration info
183:      Then set runtime options, which can override these defaults.
184:      For example,
185:           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
186:      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
187:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

189:   TSSetMaxSteps(ts,time_steps_max);
190:   TSSetMaxTime(ts,time_total_max);
191:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
192:   TSSetFromOptions(ts);

194:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195:      Solve the problem
196:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

198:   /*
199:      Evaluate initial conditions
200:   */
201:   InitialConditions(u,&appctx);

203:   /*
204:      Run the timestepping solver
205:   */
206:   TSSolve(ts,u);
207:   TSGetSolveTime(ts,&ftime);
208:   TSGetStepNumber(ts,&steps);

210:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211:      View timestepping solver info
212:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

214:   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));
215:   TSView(ts,PETSC_VIEWER_STDOUT_SELF);

217:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218:      Free work space.  All PETSc objects should be destroyed when they
219:      are no longer needed.
220:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

222:   TSDestroy(&ts);
223:   MatDestroy(&A);
224:   VecDestroy(&u);
225:   PetscViewerDestroy(&appctx.viewer1);
226:   PetscViewerDestroy(&appctx.viewer2);
227:   VecDestroy(&appctx.solution);

229:   /*
230:      Always call PetscFinalize() before exiting a program.  This routine
231:        - finalizes the PETSc libraries as well as MPI
232:        - provides summary and diagnostic information if certain runtime
233:          options are chosen (e.g., -log_view).
234:   */
235:   PetscFinalize();
236:   return 0;
237: }
238: /* --------------------------------------------------------------------- */
239: /*
240:    InitialConditions - Computes the solution at the initial time.

242:    Input Parameter:
243:    u - uninitialized solution vector (global)
244:    appctx - user-defined application context

246:    Output Parameter:
247:    u - vector with solution at initial time (global)
248: */
249: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
250: {
251:   PetscScalar    *u_localptr,h = appctx->h;
252:   PetscInt       i;

254:   /*
255:     Get a pointer to vector data.
256:     - For default PETSc vectors, VecGetArray() returns a pointer to
257:       the data array.  Otherwise, the routine is implementation dependent.
258:     - You MUST call VecRestoreArray() when you no longer need access to
259:       the array.
260:     - Note that the Fortran interface to VecGetArray() differs from the
261:       C version.  See the users manual for details.
262:   */
263:   VecGetArray(u,&u_localptr);

265:   /*
266:      We initialize the solution array by simply writing the solution
267:      directly into the array locations.  Alternatively, we could use
268:      VecSetValues() or VecSetValuesLocal().
269:   */
270:   for (i=0; i<appctx->m; i++) u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h);

272:   /*
273:      Restore vector
274:   */
275:   VecRestoreArray(u,&u_localptr);

277:   /*
278:      Print debugging information if desired
279:   */
280:   if (appctx->debug) {
281:     printf("initial guess vector\n");
282:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
283:   }

285:   return 0;
286: }
287: /* --------------------------------------------------------------------- */
288: /*
289:    ExactSolution - Computes the exact solution at a given time.

291:    Input Parameters:
292:    t - current time
293:    solution - vector in which exact solution will be computed
294:    appctx - user-defined application context

296:    Output Parameter:
297:    solution - vector with the newly computed exact solution
298: */
299: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
300: {
301:   PetscScalar    *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
302:   PetscInt       i;

304:   /*
305:      Get a pointer to vector data.
306:   */
307:   VecGetArray(solution,&s_localptr);

309:   /*
310:      Simply write the solution directly into the array locations.
311:      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
312:   */
313:   ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
314:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
315:   for (i=0; i<appctx->m; i++) s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2;

317:   /*
318:      Restore vector
319:   */
320:   VecRestoreArray(solution,&s_localptr);
321:   return 0;
322: }
323: /* --------------------------------------------------------------------- */
324: /*
325:    Monitor - User-provided routine to monitor the solution computed at
326:    each timestep.  This example plots the solution and computes the
327:    error in two different norms.

329:    Input Parameters:
330:    ts     - the timestep context
331:    step   - the count of the current step (with 0 meaning the
332:              initial condition)
333:    time   - the current time
334:    u      - the solution at this timestep
335:    ctx    - the user-provided context for this monitoring routine.
336:             In this case we use the application context which contains
337:             information about the problem size, workspace and the exact
338:             solution.
339: */
340: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
341: {
342:   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
343:   PetscReal      norm_2,norm_max;

345:   /*
346:      View a graph of the current iterate
347:   */
348:   VecView(u,appctx->viewer2);

350:   /*
351:      Compute the exact solution
352:   */
353:   ExactSolution(time,appctx->solution,appctx);

355:   /*
356:      Print debugging information if desired
357:   */
358:   if (appctx->debug) {
359:     printf("Computed solution vector\n");
360:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
361:     printf("Exact solution vector\n");
362:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
363:   }

365:   /*
366:      Compute the 2-norm and max-norm of the error
367:   */
368:   VecAXPY(appctx->solution,-1.0,u);
369:   VecNorm(appctx->solution,NORM_2,&norm_2);
370:   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
371:   VecNorm(appctx->solution,NORM_MAX,&norm_max);
372:   if (norm_2   < 1e-14) norm_2   = 0;
373:   if (norm_max < 1e-14) norm_max = 0;

375:   PetscPrintf(PETSC_COMM_WORLD,"Timestep %D: time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)time,(double)norm_2,(double)norm_max);
376:   appctx->norm_2   += norm_2;
377:   appctx->norm_max += norm_max;

379:   /*
380:      View a graph of the error
381:   */
382:   VecView(appctx->solution,appctx->viewer1);

384:   /*
385:      Print debugging information if desired
386:   */
387:   if (appctx->debug) {
388:     printf("Error vector\n");
389:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
390:   }

392:   return 0;
393: }
394: /* --------------------------------------------------------------------- */
395: /*
396:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
397:    matrix for the heat equation.

399:    Input Parameters:
400:    ts - the TS context
401:    t - current time
402:    global_in - global input vector
403:    dummy - optional user-defined context, as set by TSetRHSJacobian()

405:    Output Parameters:
406:    AA - Jacobian matrix
407:    BB - optionally different preconditioning matrix
408:    str - flag indicating matrix structure

410:   Notes:
411:   Recall that MatSetValues() uses 0-based row and column numbers
412:   in Fortran as well as in C.
413: */
414: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
415: {
416:   Mat            A       = AA;                /* Jacobian matrix */
417:   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
418:   PetscInt       mstart  = 0;
419:   PetscInt       mend    = appctx->m;
420:   PetscInt       i,idx[3];
421:   PetscScalar    v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;

423:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
424:      Compute entries for the locally owned part of the matrix
425:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
426:   /*
427:      Set matrix rows corresponding to boundary data
428:   */

430:   mstart = 0;
431:   v[0]   = 1.0;
432:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
433:   mstart++;

435:   mend--;
436:   v[0] = 1.0;
437:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

439:   /*
440:      Set matrix rows corresponding to interior data.  We construct the
441:      matrix one row at a time.
442:   */
443:   v[0] = sone; v[1] = stwo; v[2] = sone;
444:   for (i=mstart; i<mend; i++) {
445:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
446:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
447:   }

449:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
450:      Complete the matrix assembly process and set some options
451:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
452:   /*
453:      Assemble matrix, using the 2-step process:
454:        MatAssemblyBegin(), MatAssemblyEnd()
455:      Computations can be done while messages are in transition
456:      by placing code between these two statements.
457:   */
458:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
459:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

461:   /*
462:      Set and option to indicate that we will never add a new nonzero location
463:      to the matrix. If we do, it will generate an error.
464:   */
465:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

467:   return 0;
468: }

470: /*TEST

472:     test:
473:       requires: x

475:     test:
476:       suffix: nox
477:       args: -nox

479: TEST*/