Actual source code: ex5.c

petsc-3.7.5 2017-01-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:    Processors: 1
 14: */

 16: /* ------------------------------------------------------------------------

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

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

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

 38:    Notes:
 39:    This code demonstrates the TS solver interface to two variants of
 40:    linear problems, u_t = f(u,t), namely
 41:      - time-dependent f:   f(u,t) is a function of t
 42:      - time-independent f: f(u,t) is simply just f(u)

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

 46:   ------------------------------------------------------------------------- */

 48: /*
 49:    Include "petscts.h" so that we can use TS solvers.  Note that this file
 50:    automatically includes:
 51:      petscsys.h       - base PETSc routines   petscvec.h  - vectors
 52:      petscmat.h  - matrices
 53:      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
 54:      petscviewer.h - viewers               petscpc.h   - preconditioners
 55:      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
 56: */
 57: #include <petscts.h>
 58: #include <petscdraw.h>

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

 73: /*
 74:    User-defined routines
 75: */
 76: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
 77: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
 78: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 79: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);

 83: int main(int argc,char **argv)
 84: {
 85:   AppCtx         appctx;                 /* user-defined application context */
 86:   TS             ts;                     /* timestepping context */
 87:   Mat            A;                      /* matrix data structure */
 88:   Vec            u;                      /* approximate solution vector */
 89:   PetscReal      time_total_max = 100.0; /* default max total time */
 90:   PetscInt       time_steps_max = 100;   /* default max timesteps */
 91:   PetscDraw      draw;                   /* drawing context */
 93:   PetscInt       steps,m;
 94:   PetscMPIInt    size;
 95:   PetscBool      flg;
 96:   PetscReal      dt,ftime;

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 99:      Initialize program and set problem parameters
100:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

106:   m               = 60;
107:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
108:   PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);
109:   appctx.m        = m;
110:   appctx.h        = 1.0/(m-1.0);
111:   appctx.norm_2   = 0.0;
112:   appctx.norm_max = 0.0;

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

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:      Create vector data structures
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

120:   /*
121:      Create vector data structures for approximate and exact solutions
122:   */
123:   VecCreateSeq(PETSC_COMM_SELF,m,&u);
124:   VecDuplicate(u,&appctx.solution);

126:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127:      Set up displays to show graphs of the solution and error
128:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

130:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
131:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
132:   PetscDrawSetDoubleBuffer(draw);
133:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
134:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
135:   PetscDrawSetDoubleBuffer(draw);

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Create timestepping solver context
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

141:   TSCreate(PETSC_COMM_SELF,&ts);
142:   TSSetProblemType(ts,TS_LINEAR);

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:      Set optional user-defined monitoring routine
146:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

152:      Create matrix data structure; set matrix evaluation routine.
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

155:   MatCreate(PETSC_COMM_SELF,&A);
156:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
157:   MatSetFromOptions(A);
158:   MatSetUp(A);

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

181:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182:      Set solution vector and initial timestep
183:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

185:   dt   = appctx.h*appctx.h/2.0;
186:   TSSetInitialTimeStep(ts,0.0,dt);
187:   TSSetSolution(ts,u);

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

199:   TSSetDuration(ts,time_steps_max,time_total_max);
200:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
201:   TSSetFromOptions(ts);

203:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204:      Solve the problem
205:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

207:   /*
208:      Evaluate initial conditions
209:   */
210:   InitialConditions(u,&appctx);

212:   /*
213:      Run the timestepping solver
214:   */
215:   TSSolve(ts,u);
216:   TSGetSolveTime(ts,&ftime);
217:   TSGetTimeStepNumber(ts,&steps);

219:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220:      View timestepping solver info
221:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

226:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227:      Free work space.  All PETSc objects should be destroyed when they
228:      are no longer needed.
229:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

231:   TSDestroy(&ts);
232:   MatDestroy(&A);
233:   VecDestroy(&u);
234:   PetscViewerDestroy(&appctx.viewer1);
235:   PetscViewerDestroy(&appctx.viewer2);
236:   VecDestroy(&appctx.solution);

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

253:    Input Parameter:
254:    u - uninitialized solution vector (global)
255:    appctx - user-defined application context

257:    Output Parameter:
258:    u - vector with solution at initial time (global)
259: */
260: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
261: {
262:   PetscScalar    *u_localptr,h = appctx->h;
263:   PetscInt       i;

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

277:   /*
278:      We initialize the solution array by simply writing the solution
279:      directly into the array locations.  Alternatively, we could use
280:      VecSetValues() or VecSetValuesLocal().
281:   */
282:   for (i=0; i<appctx->m; i++) u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h);

284:   /*
285:      Restore vector
286:   */
287:   VecRestoreArray(u,&u_localptr);

289:   /*
290:      Print debugging information if desired
291:   */
292:   if (appctx->debug) {
293:     printf("initial guess vector\n");
294:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
295:   }

297:   return 0;
298: }
299: /* --------------------------------------------------------------------- */
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,tc = t;
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 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
329:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
330:   for (i=0; i<appctx->m; i++) s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2;

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

346:    Input Parameters:
347:    ts     - the timestep context
348:    step   - the count of the current step (with 0 meaning the
349:              initial condition)
350:    time   - 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 time,Vec u,void *ctx)
358: {
359:   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
361:   PetscReal      norm_2,norm_max;

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

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

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

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

393:   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);
394:   appctx->norm_2   += norm_2;
395:   appctx->norm_max += norm_max;

397:   /*
398:      View a graph of the error
399:   */
400:   VecView(appctx->solution,appctx->viewer1);

402:   /*
403:      Print debugging information if desired
404:   */
405:   if (appctx->debug) {
406:     printf("Error vector\n");
407:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
408:   }

410:   return 0;
411: }
412: /* --------------------------------------------------------------------- */
415: /*
416:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
417:    matrix for the heat equation.

419:    Input Parameters:
420:    ts - the TS context
421:    t - current time
422:    global_in - global input vector
423:    dummy - optional user-defined context, as set by TSetRHSJacobian()

425:    Output Parameters:
426:    AA - Jacobian matrix
427:    BB - optionally different preconditioning matrix
428:    str - flag indicating matrix structure

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

444:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
445:      Compute entries for the locally owned part of the matrix
446:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
447:   /*
448:      Set matrix rows corresponding to boundary data
449:   */

451:   mstart = 0;
452:   v[0]   = 1.0;
453:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
454:   mstart++;

456:   mend--;
457:   v[0] = 1.0;
458:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

460:   /*
461:      Set matrix rows corresponding to interior data.  We construct the
462:      matrix one row at a time.
463:   */
464:   v[0] = sone; v[1] = stwo; v[2] = sone;
465:   for (i=mstart; i<mend; i++) {
466:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
467:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
468:   }

470:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
471:      Complete the matrix assembly process and set some options
472:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
473:   /*
474:      Assemble matrix, using the 2-step process:
475:        MatAssemblyBegin(), MatAssemblyEnd()
476:      Computations can be done while messages are in transition
477:      by placing code between these two statements.
478:   */
479:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
480:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

482:   /*
483:      Set and option to indicate that we will never add a new nonzero location
484:      to the matrix. If we do, it will generate an error.
485:   */
486:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

488:   return 0;
489: }