Actual source code: ex6.c

petsc-3.7.3 2016-08-01
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,NULL,"-m",&m,NULL);
114:   PetscOptionsHasName(NULL,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,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:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
208:   TSSetFromOptions(ts);

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

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

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

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

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

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

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

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

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

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

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

284:   /*
285:      We initialize the solution array by simply writing the solution
286:      directly into the array locations.  Alternatively, we could use
287:      VecSetValues() or VecSetValuesLocal().
288:   */
289:   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);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

400:   TSGetTimeStep(ts,&dt);
401:   if (norm_2 > 1.e-2) {
402:     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);
403:   }
404:   appctx->norm_2   += norm_2;
405:   appctx->norm_max += norm_max;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

530:   return 0;
531: }