Actual source code: ex6.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  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:    Routines: TSCreate(); TSSetSolution(); TSSetRHSJacobian(), TSSetIJacobian();
 14:    Routines: TSSetInitialTimeStep(); TSSetDuration(); TSMonitorSet();
 15:    Routines: TSSetFromOptions(); TSStep(); TSDestroy();
 16:    Routines: TSSetTimeStep(); TSGetTimeStep();
 17:    Processors: 1
 18: */

 20: /* ------------------------------------------------------------------------

 22:    This program solves the one-dimensional heat equation (also called the
 23:    diffusion equation),
 24:        u_t = u_xx,
 25:    on the domain 0 <= x <= 1, with the boundary conditions
 26:        u(t,0) = 0, u(t,1) = 0,
 27:    and the initial condition
 28:        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
 29:    This is a linear, second-order, parabolic equation.

 31:    We discretize the right-hand side using finite differences with
 32:    uniform grid spacing h:
 33:        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
 34:    We then demonstrate time evolution using the various TS methods by
 35:    running the program via
 36:        ex3 -ts_type <timestepping solver>

 38:    We compare the approximate solution with the exact solution, given by
 39:        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
 40:                       3*exp(-4*pi*pi*t) * sin(2*pi*x)

 42:    Notes:
 43:    This code demonstrates the TS solver interface to two variants of
 44:    linear problems, u_t = f(u,t), namely
 45:      - time-dependent f:   f(u,t) is a function of t
 46:      - time-independent f: f(u,t) is simply f(u)

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

 50:   ------------------------------------------------------------------------- */

 52: /*
 53:    Include "ts.h" so that we can use TS solvers.  Note that this file
 54:    automatically includes:
 55:      petscsys.h  - base PETSc routines   vec.h  - vectors
 56:      sys.h    - system routines       mat.h  - matrices
 57:      is.h     - index sets            ksp.h  - Krylov subspace methods
 58:      viewer.h - viewers               pc.h   - preconditioners
 59:      snes.h - nonlinear solvers
 60: */

 62: #include <petscts.h>
 63: #include <petscdraw.h>

 65: /*
 66:    User-defined application context - contains data needed by the
 67:    application-provided call-back routines.
 68: */
 69: typedef struct {
 70:   Vec         solution;          /* global exact solution vector */
 71:   PetscInt    m;                 /* total number of grid points */
 72:   PetscReal   h;                 /* mesh width h = 1/(m-1) */
 73:   PetscBool   debug;             /* flag (1 indicates activation of debugging printouts) */
 74:   PetscViewer viewer1, viewer2;  /* viewers for the solution and error */
 75:   PetscReal   norm_2, norm_max;  /* error norms */
 76: } AppCtx;

 78: /*
 79:    User-defined routines
 80: */
 81: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
 82: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
 83: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 84: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
 85: extern PetscErrorCode MyBCRoutine(TS,PetscReal,Vec,void*);

 89: int main(int argc,char **argv)
 90: {
 91:   AppCtx         appctx;                 /* user-defined application context */
 92:   TS             ts;                     /* timestepping context */
 93:   Mat            A;                      /* matrix data structure */
 94:   Vec            u;                      /* approximate solution vector */
 95:   PetscReal      time_total_max = 100.0; /* default max total time */
 96:   PetscInt       time_steps_max = 100;   /* default max timesteps */
 97:   PetscDraw      draw;                   /* drawing context */
 99:   PetscInt       steps, m;
100:   PetscMPIInt    size;
101:   PetscReal      dt;
102:   PetscReal      ftime;
103:   PetscBool      flg;
104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:      Initialize program and set problem parameters
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   PetscInitialize(&argc,&argv,(char*)0,help);
109:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
110:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

112:   m    = 60;
113:   PetscOptionsGetInt(NULL,"-m",&m,NULL);
114:   PetscOptionsHasName(NULL,"-debug",&appctx.debug);

116:   appctx.m        = m;
117:   appctx.h        = 1.0/(m-1.0);
118:   appctx.norm_2   = 0.0;
119:   appctx.norm_max = 0.0;

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

123:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124:      Create vector data structures
125:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

127:   /*
128:      Create vector data structures for approximate and exact solutions
129:   */
130:   VecCreateSeq(PETSC_COMM_SELF,m,&u);
131:   VecDuplicate(u,&appctx.solution);

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:      Set up displays to show graphs of the solution and error
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

137:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
138:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
139:   PetscDrawSetDoubleBuffer(draw);
140:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
141:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
142:   PetscDrawSetDoubleBuffer(draw);

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:      Create timestepping solver context
146:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

148:   TSCreate(PETSC_COMM_SELF,&ts);
149:   TSSetProblemType(ts,TS_LINEAR);

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:      Set optional user-defined monitoring routine
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

157:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

159:      Create matrix data structure; set matrix evaluation routine.
160:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

162:   MatCreate(PETSC_COMM_SELF,&A);
163:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
164:   MatSetFromOptions(A);
165:   MatSetUp(A);

167:   PetscOptionsHasName(NULL,"-time_dependent_rhs",&flg);
168:   if (flg) {
169:     /*
170:        For linear problems with a time-dependent f(u,t) in the equation
171:        u_t = f(u,t), the user provides the discretized right-hand-side
172:        as a time-dependent matrix.
173:     */
174:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
175:     TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
176:   } else {
177:     /*
178:        For linear problems with a time-independent f(u) in the equation
179:        u_t = f(u), the user provides the discretized right-hand-side
180:        as a matrix only once, and then sets a null matrix evaluation
181:        routine.
182:     */
183:     RHSMatrixHeat(ts,0.0,u,A,A,&appctx);
184:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
185:     TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
186:   }

188:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189:      Set solution vector and initial timestep
190:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

192:   dt   = appctx.h*appctx.h/2.0;
193:   TSSetInitialTimeStep(ts,0.0,dt);
194:   TSSetSolution(ts,u);

196:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197:      Customize timestepping solver:
198:        - Set the solution method to be the Backward Euler method.
199:        - Set timestepping duration info
200:      Then set runtime options, which can override these defaults.
201:      For example,
202:           -ts_max_steps <maxsteps> -ts_final_time <maxtime>
203:      to override the defaults set by TSSetDuration().
204:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

206:   TSSetDuration(ts,time_steps_max,time_total_max);
207:   TSSetFromOptions(ts);

209:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210:      Solve the problem
211:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

213:   /*
214:      Evaluate initial conditions
215:   */
216:   InitialConditions(u,&appctx);

218:   /*
219:      Run the timestepping solver
220:   */
221:   TSSolve(ts,u);
222:   TSGetSolveTime(ts,&ftime);
223:   TSGetTimeStepNumber(ts,&steps);

225:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226:      View timestepping solver info
227:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

232:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233:      Free work space.  All PETSc objects should be destroyed when they
234:      are no longer needed.
235:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

237:   TSDestroy(&ts);
238:   MatDestroy(&A);
239:   VecDestroy(&u);
240:   PetscViewerDestroy(&appctx.viewer1);
241:   PetscViewerDestroy(&appctx.viewer2);
242:   VecDestroy(&appctx.solution);

244:   /*
245:      Always call PetscFinalize() before exiting a program.  This routine
246:        - finalizes the PETSc libraries as well as MPI
247:        - provides summary and diagnostic information if certain runtime
248:          options are chosen (e.g., -log_summary).
249:   */
250:   PetscFinalize();
251:   return 0;
252: }
253: /* --------------------------------------------------------------------- */
256: /*
257:    InitialConditions - Computes the solution at the initial time.

259:    Input Parameter:
260:    u - uninitialized solution vector (global)
261:    appctx - user-defined application context

263:    Output Parameter:
264:    u - vector with solution at initial time (global)
265: */
266: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
267: {
268:   PetscScalar    *u_localptr;
269:   PetscInt       i;

272:   /*
273:     Get a pointer to vector data.
274:     - For default PETSc vectors, VecGetArray() returns a pointer to
275:       the data array.  Otherwise, the routine is implementation dependent.
276:     - You MUST call VecRestoreArray() when you no longer need access to
277:       the array.
278:     - Note that the Fortran interface to VecGetArray() differs from the
279:       C version.  See the users manual for details.
280:   */
281:   VecGetArray(u,&u_localptr);

283:   /*
284:      We initialize the solution array by simply writing the solution
285:      directly into the array locations.  Alternatively, we could use
286:      VecSetValues() or VecSetValuesLocal().
287:   */
288:   for (i=0; i<appctx->m; i++) u_localptr[i] = PetscSinReal(PETSC_PI*i*6.*appctx->h) + 3.*PetscSinReal(PETSC_PI*i*2.*appctx->h);

290:   /*
291:      Restore vector
292:   */
293:   VecRestoreArray(u,&u_localptr);

295:   /*
296:      Print debugging information if desired
297:   */
298:   if (appctx->debug) {
299:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
300:   }

302:   return 0;
303: }
304: /* --------------------------------------------------------------------- */
307: /*
308:    ExactSolution - Computes the exact solution at a given time.

310:    Input Parameters:
311:    t - current time
312:    solution - vector in which exact solution will be computed
313:    appctx - user-defined application context

315:    Output Parameter:
316:    solution - vector with the newly computed exact solution
317: */
318: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
319: {
320:   PetscScalar    *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
321:   PetscInt       i;

324:   /*
325:      Get a pointer to vector data.
326:   */
327:   VecGetArray(solution,&s_localptr);

329:   /*
330:      Simply write the solution directly into the array locations.
331:      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
332:   */
333:   ex1 = PetscExpReal(-36.*PETSC_PI*PETSC_PI*t); ex2 = PetscExpReal(-4.*PETSC_PI*PETSC_PI*t);
334:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
335:   for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinReal(PetscRealPart(sc1)*(PetscReal)i)*ex1 + 3.*PetscSinReal(PetscRealPart(sc2)*(PetscReal)i)*ex2;

337:   /*
338:      Restore vector
339:   */
340:   VecRestoreArray(solution,&s_localptr);
341:   return 0;
342: }
343: /* --------------------------------------------------------------------- */
346: /*
347:    Monitor - User-provided routine to monitor the solution computed at
348:    each timestep.  This example plots the solution and computes the
349:    error in two different norms.

351:    This example also demonstrates changing the timestep via TSSetTimeStep().

353:    Input Parameters:
354:    ts     - the timestep context
355:    step   - the count of the current step (with 0 meaning the
356:              initial condition)
357:    crtime  - the current time
358:    u      - the solution at this timestep
359:    ctx    - the user-provided context for this monitoring routine.
360:             In this case we use the application context which contains
361:             information about the problem size, workspace and the exact
362:             solution.
363: */
364: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal crtime,Vec u,void *ctx)
365: {
366:   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
368:   PetscReal      norm_2, norm_max, dt, dttol;
369:   PetscBool      flg;

371:   /*
372:      View a graph of the current iterate
373:   */
374:   VecView(u,appctx->viewer2);

376:   /*
377:      Compute the exact solution
378:   */
379:   ExactSolution(crtime,appctx->solution,appctx);

381:   /*
382:      Print debugging information if desired
383:   */
384:   if (appctx->debug) {
385:     PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n");
386:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
387:     PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");
388:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
389:   }

391:   /*
392:      Compute the 2-norm and max-norm of the error
393:   */
394:   VecAXPY(appctx->solution,-1.0,u);
395:   VecNorm(appctx->solution,NORM_2,&norm_2);
396:   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
397:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

399:   TSGetTimeStep(ts,&dt);
400:   if (norm_2 > 1.e-2) {
401:     PetscPrintf(PETSC_COMM_SELF,"Timestep %D: step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)dt,(double)crtime,(double)norm_2,(double)norm_max);
402:   }
403:   appctx->norm_2   += norm_2;
404:   appctx->norm_max += norm_max;

406:   dttol = .0001;
407:   PetscOptionsGetReal(NULL,"-dttol",&dttol,&flg);
408:   if (dt < dttol) {
409:     dt  *= .999;
410:     TSSetTimeStep(ts,dt);
411:   }

413:   /*
414:      View a graph of the error
415:   */
416:   VecView(appctx->solution,appctx->viewer1);

418:   /*
419:      Print debugging information if desired
420:   */
421:   if (appctx->debug) {
422:     PetscPrintf(PETSC_COMM_SELF,"Error vector\n");
423:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
424:   }

426:   return 0;
427: }
428: /* --------------------------------------------------------------------- */
431: /*
432:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
433:    matrix for the heat equation.

435:    Input Parameters:
436:    ts - the TS context
437:    t - current time
438:    global_in - global input vector
439:    dummy - optional user-defined context, as set by TSetRHSJacobian()

441:    Output Parameters:
442:    AA - Jacobian matrix
443:    BB - optionally different preconditioning matrix
444:    str - flag indicating matrix structure

446:    Notes:
447:    Recall that MatSetValues() uses 0-based row and column numbers
448:    in Fortran as well as in C.
449: */
450: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
451: {
452:   Mat            A       = AA;                /* Jacobian matrix */
453:   AppCtx         *appctx = (AppCtx*) ctx;      /* user-defined application context */
454:   PetscInt       mstart  = 0;
455:   PetscInt       mend    = appctx->m;
457:   PetscInt       i, idx[3];
458:   PetscScalar    v[3], stwo = -2./(appctx->h*appctx->h), sone = -.5*stwo;

460:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
461:      Compute entries for the locally owned part of the matrix
462:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
463:   /*
464:      Set matrix rows corresponding to boundary data
465:   */

467:   mstart = 0;
468:   v[0]   = 1.0;
469:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
470:   mstart++;

472:   mend--;
473:   v[0] = 1.0;
474:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

476:   /*
477:      Set matrix rows corresponding to interior data.  We construct the
478:      matrix one row at a time.
479:   */
480:   v[0] = sone; v[1] = stwo; v[2] = sone;
481:   for (i=mstart; i<mend; i++) {
482:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
483:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
484:   }

486:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
487:      Complete the matrix assembly process and set some options
488:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
489:   /*
490:      Assemble matrix, using the 2-step process:
491:        MatAssemblyBegin(), MatAssemblyEnd()
492:      Computations can be done while messages are in transition
493:      by placing code between these two statements.
494:   */
495:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
496:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

498:   /*
499:      Set and option to indicate that we will never add a new nonzero location
500:      to the matrix. If we do, it will generate an error.
501:   */
502:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

504:   return 0;
505: }
506: /* --------------------------------------------------------------------- */
509: /*
510:    Input Parameters:
511:    ts - the TS context
512:    t - current time
513:    f - function
514:    ctx - optional user-defined context, as set by TSetBCFunction()
515:  */
516: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
517: {
518:   AppCtx         *appctx = (AppCtx*) ctx;      /* user-defined application context */
520:   PetscInt       m = appctx->m;
521:   PetscScalar    *fa;

523:   VecGetArray(f,&fa);
524:   fa[0]   = 0.0;
525:   fa[m-1] = 1.0;
526:   VecRestoreArray(f,&fa);
527:   PetscPrintf(PETSC_COMM_SELF,"t=%g\n",(double)t);

529:   return 0;
530: }