Actual source code: ex6.c

petsc-master 2017-04-20
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*);

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

186:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187:      Set solution vector and initial timestep
188:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

204:   TSSetDuration(ts,time_steps_max,time_total_max);
205:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
206:   TSSetFromOptions(ts);

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

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

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

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

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

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

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

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

256:    Input Parameter:
257:    u - uninitialized solution vector (global)
258:    appctx - user-defined application context

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

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

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

287:   /*
288:      Restore vector
289:   */
290:   VecRestoreArray(u,&u_localptr);

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

299:   return 0;
300: }
301: /* --------------------------------------------------------------------- */
302: /*
303:    ExactSolution - Computes the exact solution at a given time.

305:    Input Parameters:
306:    t - current time
307:    solution - vector in which exact solution will be computed
308:    appctx - user-defined application context

310:    Output Parameter:
311:    solution - vector with the newly computed exact solution
312: */
313: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
314: {
315:   PetscScalar    *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
316:   PetscInt       i;

319:   /*
320:      Get a pointer to vector data.
321:   */
322:   VecGetArray(solution,&s_localptr);

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

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

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

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

364:   /*
365:      View a graph of the current iterate
366:   */
367:   VecView(u,appctx->viewer2);

369:   /*
370:      Compute the exact solution
371:   */
372:   ExactSolution(crtime,appctx->solution,appctx);

374:   /*
375:      Print debugging information if desired
376:   */
377:   if (appctx->debug) {
378:     PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n");
379:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
380:     PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");
381:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
382:   }

384:   /*
385:      Compute the 2-norm and max-norm of the error
386:   */
387:   VecAXPY(appctx->solution,-1.0,u);
388:   VecNorm(appctx->solution,NORM_2,&norm_2);
389:   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
390:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

392:   TSGetTimeStep(ts,&dt);
393:   if (norm_2 > 1.e-2) {
394:     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);
395:   }
396:   appctx->norm_2   += norm_2;
397:   appctx->norm_max += norm_max;

399:   dttol = .0001;
400:   PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,&flg);
401:   if (dt < dttol) {
402:     dt  *= .999;
403:     TSSetTimeStep(ts,dt);
404:   }

406:   /*
407:      View a graph of the error
408:   */
409:   VecView(appctx->solution,appctx->viewer1);

411:   /*
412:      Print debugging information if desired
413:   */
414:   if (appctx->debug) {
415:     PetscPrintf(PETSC_COMM_SELF,"Error vector\n");
416:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
417:   }

419:   return 0;
420: }
421: /* --------------------------------------------------------------------- */
422: /*
423:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
424:    matrix for the heat equation.

426:    Input Parameters:
427:    ts - the TS context
428:    t - current time
429:    global_in - global input vector
430:    dummy - optional user-defined context, as set by TSetRHSJacobian()

432:    Output Parameters:
433:    AA - Jacobian matrix
434:    BB - optionally different preconditioning matrix
435:    str - flag indicating matrix structure

437:    Notes:
438:    Recall that MatSetValues() uses 0-based row and column numbers
439:    in Fortran as well as in C.
440: */
441: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
442: {
443:   Mat            A       = AA;                /* Jacobian matrix */
444:   AppCtx         *appctx = (AppCtx*) ctx;      /* user-defined application context */
445:   PetscInt       mstart  = 0;
446:   PetscInt       mend    = appctx->m;
448:   PetscInt       i, idx[3];
449:   PetscScalar    v[3], stwo = -2./(appctx->h*appctx->h), sone = -.5*stwo;

451:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
452:      Compute entries for the locally owned part of the matrix
453:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
454:   /*
455:      Set matrix rows corresponding to boundary data
456:   */

458:   mstart = 0;
459:   v[0]   = 1.0;
460:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
461:   mstart++;

463:   mend--;
464:   v[0] = 1.0;
465:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

467:   /*
468:      Set matrix rows corresponding to interior data.  We construct the
469:      matrix one row at a time.
470:   */
471:   v[0] = sone; v[1] = stwo; v[2] = sone;
472:   for (i=mstart; i<mend; i++) {
473:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
474:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
475:   }

477:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
478:      Complete the matrix assembly process and set some options
479:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
480:   /*
481:      Assemble matrix, using the 2-step process:
482:        MatAssemblyBegin(), MatAssemblyEnd()
483:      Computations can be done while messages are in transition
484:      by placing code between these two statements.
485:   */
486:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
487:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

489:   /*
490:      Set and option to indicate that we will never add a new nonzero location
491:      to the matrix. If we do, it will generate an error.
492:   */
493:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

495:   return 0;
496: }
497: /* --------------------------------------------------------------------- */
498: /*
499:    Input Parameters:
500:    ts - the TS context
501:    t - current time
502:    f - function
503:    ctx - optional user-defined context, as set by TSetBCFunction()
504:  */
505: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
506: {
507:   AppCtx         *appctx = (AppCtx*) ctx;      /* user-defined application context */
509:   PetscInt       m = appctx->m;
510:   PetscScalar    *fa;

512:   VecGetArray(f,&fa);
513:   fa[0]   = 0.0;
514:   fa[m-1] = 1.0;
515:   VecRestoreArray(f,&fa);
516:   PetscPrintf(PETSC_COMM_SELF,"t=%g\n",(double)t);

518:   return 0;
519: }