Actual source code: ex6.c

petsc-3.4.5 2014-06-29
  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*,MatStructure*,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:     MatStructure A_structure;
184:     RHSMatrixHeat(ts,0.0,u,&A,&A,&A_structure,&appctx);
185:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
186:     TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
187:   }

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

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

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

207:   TSSetDuration(ts,time_steps_max,time_total_max);
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",
231:               appctx.norm_2/steps,appctx.norm_max/steps);
232:   TSView(ts,PETSC_VIEWER_STDOUT_SELF);

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

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

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

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

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

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

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

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

297:   /*
298:      Print debugging information if desired
299:   */
300:   if (appctx->debug) {
301:      printf("initial guess vector\n");
302:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
303:   }

305:   return 0;
306: }
307: /* --------------------------------------------------------------------- */
310: /*
311:    ExactSolution - Computes the exact solution at a given time.

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

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

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

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

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

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

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

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

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

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

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

402:   TSGetTimeStep(ts,&dt);
403:   if (norm_2 > 1.e-2) {
404:     printf("Timestep %d: step size = %G, time = %G, 2-norm error = %G, max norm error = %G\n",
405:          (int)step,dt,crtime,norm_2,norm_max);
406:   }
407:   appctx->norm_2   += norm_2;
408:   appctx->norm_max += norm_max;

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

417:   /*
418:      View a graph of the error
419:   */
420:   VecView(appctx->solution,appctx->viewer1);

422:   /*
423:      Print debugging information if desired
424:   */
425:   if (appctx->debug) {
426:     printf("Error vector\n");
427:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
428:   }

430:   return 0;
431: }
432: /* --------------------------------------------------------------------- */
435: /*
436:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
437:    matrix for the heat equation.

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

445:    Output Parameters:
446:    AA - Jacobian matrix
447:    BB - optionally different preconditioning matrix
448:    str - flag indicating matrix structure

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

464:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
465:      Compute entries for the locally owned part of the matrix
466:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
467:   /*
468:      Set matrix rows corresponding to boundary data
469:   */

471:   mstart = 0;
472:   v[0]   = 1.0;
473:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
474:   mstart++;

476:   mend--;
477:   v[0] = 1.0;
478:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

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

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

502:   /*
503:      Set flag to indicate that the Jacobian matrix retains an identical
504:      nonzero structure throughout all timestepping iterations (although the
505:      values of the entries change). Thus, we can save some work in setting
506:      up the preconditioner (e.g., no need to redo symbolic factorization for
507:      ILU/ICC preconditioners).
508:       - If the nonzero structure of the matrix is different during
509:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
510:         must be used instead.  If you are unsure whether the matrix
511:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
512:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
513:         believes your assertion and does not check the structure
514:         of the matrix.  If you erroneously claim that the structure
515:         is the same when it actually is not, the new preconditioner
516:         will not function correctly.  Thus, use this optimization
517:         feature with caution!
518:   */
519:   *str = SAME_NONZERO_PATTERN;

521:   /*
522:      Set and option to indicate that we will never add a new nonzero location
523:      to the matrix. If we do, it will generate an error.
524:   */
525:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

527:   return 0;
528: }
529: /* --------------------------------------------------------------------- */
532: /*
533:    Input Parameters:
534:    ts - the TS context
535:    t - current time
536:    f - function
537:    ctx - optional user-defined context, as set by TSetBCFunction()
538:  */
539: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
540: {
541:   AppCtx         *appctx = (AppCtx*) ctx;      /* user-defined application context */
543:   PetscInt       m = appctx->m;
544:   PetscScalar    *fa;

546:   VecGetArray(f,&fa);
547:   fa[0]   = 0.0;
548:   fa[m-1] = 1.0;
549:   VecRestoreArray(f,&fa);
550:   printf("t=%g\n",t);

552:   return 0;
553: }