Actual source code: ex5.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:    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,"-m",&m,NULL);
108:   PetscOptionsHasName(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,"-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:   TSSetFromOptions(ts);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

296:   return 0;
297: }
298: /* --------------------------------------------------------------------- */
301: /*
302:    ExactSolution - Computes the exact solution at a given time.

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

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

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

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

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

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

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

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

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

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

390:   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);
391:   appctx->norm_2   += norm_2;
392:   appctx->norm_max += norm_max;

394:   /*
395:      View a graph of the error
396:   */
397:   VecView(appctx->solution,appctx->viewer1);

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

407:   return 0;
408: }
409: /* --------------------------------------------------------------------- */
412: /*
413:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
414:    matrix for the heat equation.

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

422:    Output Parameters:
423:    AA - Jacobian matrix
424:    BB - optionally different preconditioning matrix
425:    str - flag indicating matrix structure

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

441:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
442:      Compute entries for the locally owned part of the matrix
443:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
444:   /*
445:      Set matrix rows corresponding to boundary data
446:   */

448:   mstart = 0;
449:   v[0]   = 1.0;
450:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
451:   mstart++;

453:   mend--;
454:   v[0] = 1.0;
455:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

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

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

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

485:   return 0;
486: }