Actual source code: ex3.c

  2: /* Program usage:  ex3 [-help] [all PETSc options] */

  4: static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\
  5: Input parameters include:\n\
  6:   -m <points>, where <points> = number of grid points\n\
  7:   -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
  8:   -time_dependent_bc  : Treat the problem as having time-dependent boundary conditions\n\
  9:   -debug              : Activate debugging printouts\n\
 10:   -nox                : Deactivate x-window graphics\n\n";

 12: /*
 13:    Concepts: TS^time-dependent linear problems
 14:    Concepts: TS^heat equation
 15:    Concepts: TS^diffusion equation
 16:    Processors: 1
 17: */

 19: /* ------------------------------------------------------------------------

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

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

 37:    We compare the approximate solution with the exact solution, given by
 38:        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
 39:                       3*exp(-4*pi*pi*t) * sin(2*pi*x)

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

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

 49:   ------------------------------------------------------------------------- */

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

 61:  #include petscts.h

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

 76: /* 
 77:    User-defined routines
 78: */

 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,ftime;
100:   PetscTruth     flg;

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:      Initialize program and set problem parameters
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: 
106:   PetscInitialize(&argc,&argv,(char*)0,help);
107:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
108:   if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");

110:   m    = 60;
111:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
112:   PetscOptionsHasName(PETSC_NULL,"-debug",&appctx.debug);
113:   appctx.m        = m;
114:   appctx.h        = 1.0/(m-1.0);
115:   appctx.norm_2   = 0.0;
116:   appctx.norm_max = 0.0;
117:   PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");

119:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120:      Create vector data structures
121:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130:      Set up displays to show graphs of the solution and error 
131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:      Create timestepping solver context
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

144:   TSCreate(PETSC_COMM_SELF,&ts);
145:   TSSetProblemType(ts,TS_LINEAR);

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:      Set optional user-defined monitoring routine
149:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

151:   TSSetMonitor(ts,Monitor,&appctx,PETSC_NULL);

153:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

155:      Create matrix data structure; set matrix evaluation routine.
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

158:   MatCreate(PETSC_COMM_SELF,&A);
159:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
160:   MatSetFromOptions(A);

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

182:   /* Treat the problem as having time-dependent boundary conditions */
183:   PetscOptionsHasName(PETSC_NULL,"-time_dependent_bc",&flg);
184:   if (flg) {
185:     TSSetRHSBoundaryConditions(ts,MyBCRoutine,&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_max_time <maxtime>
203:      to override the defaults set by TSSetDuration().
204:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

206:   TSSetDuration(ts,time_steps_max,time_total_max);
207:   TSSetFromOptions(ts);

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

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

218:   /*
219:      Run the timestepping solver
220:   */
221:   TSStep(ts,&steps,&ftime);

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

227:   PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %g, avg. error (max norm) = %g\n",
228:               appctx.norm_2/steps,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_summary). 
248:   */
249:   PetscFinalize();
250:   return 0;
251: }
252: /* --------------------------------------------------------------------- */
255: /*
256:    InitialConditions - Computes the solution at the initial time. 

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

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

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

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

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

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

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

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

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

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

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

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

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

358:    Input Parameters:
359:    ts     - the timestep context
360:    step   - the count of the current step (with 0 meaning the
361:              initial condition)
362:    time   - the current time
363:    u      - the solution at this timestep
364:    ctx    - the user-provided context for this monitoring routine.
365:             In this case we use the application context which contains 
366:             information about the problem size, workspace and the exact 
367:             solution.
368: */
369: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
370: {
371:   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
373:   PetscReal      norm_2,norm_max,dt,dttol;
374:   PetscScalar    mone = -1.0;
375:   /* 
376:      View a graph of the current iterate
377:   */
378:   VecView(u,appctx->viewer2);

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

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

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

403:   TSGetTimeStep(ts,&dt);
404:   PetscPrintf(PETSC_COMM_WORLD,"Timestep %3D: step size = %-11g, time = %-11g, 2-norm error = %-11g, max norm error = %-11g\n",
405:          step,dt,time,norm_2,norm_max);
406:   appctx->norm_2   += norm_2;
407:   appctx->norm_max += norm_max;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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