Actual source code: ex6.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: ex3.c,v 1.7 1999/07/21 14:35:38 curfman Exp curfman $";
  3: #endif

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

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

 15: /*
 16:    Concepts: TS^time-dependent linear problems
 17:    Concepts: TS^heat equation
 18:    Concepts: TS^diffusion equation
 19:    Routines: TSCreate(); TSSetSolution(); TSSetRHSMatrix();
 20:    Routines: TSSetInitialTimeStep(); TSSetDuration(); TSSetMonitor();
 21:    Routines: TSSetFromOptions(); TSStep(); TSDestroy(); 
 22:    Routines: TSSetTimeStep(); TSGetTimeStep();
 23:    Processors: 1
 24: */

 26: /* ------------------------------------------------------------------------

 28:    This program solves the one-dimensional heat equation (also called the
 29:    diffusion equation),
 30:        u_t = u_xx, 
 31:    on the domain 0 <= x <= 1, with the boundary conditions
 32:        u(t,0) = 0, u(t,1) = 0,
 33:    and the initial condition
 34:        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
 35:    This is a linear, second-order, parabolic equation.

 37:    We discretize the right-hand side using finite differences with
 38:    uniform grid spacing h:
 39:        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
 40:    We then demonstrate time evolution using the various TS methods by
 41:    running the program via
 42:        ex3 -ts_type <timestepping solver>

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

 48:    Notes:
 49:    This code demonstrates the TS solver interface to two variants of 
 50:    linear problems, u_t = f(u,t), namely
 51:      - time-dependent f:   f(u,t) is a function of t
 52:      - time-independent f: f(u,t) is simply f(u)

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

 56:   ------------------------------------------------------------------------- */

 58: /* 
 59:    Include "ts.h" so that we can use TS solvers.  Note that this file
 60:    automatically includes:
 61:      petsc.h  - base PETSc routines   vec.h  - vectors
 62:      sys.h    - system routines       mat.h  - matrices
 63:      is.h     - index sets            ksp.h  - Krylov subspace methods
 64:      viewer.h - viewers               pc.h   - preconditioners
 65:      snes.h - nonlinear solvers
 66: */

 68:  #include petscts.h

 70: /* 
 71:    User-defined application context - contains data needed by the 
 72:    application-provided call-back routines.
 73: */
 74: typedef struct {
 75:   Vec      solution;          /* global exact solution vector */
 76:   int      m;                 /* total number of grid points */
 77:   double   h;                 /* mesh width h = 1/(m-1) */
 78:   int      debug;             /* flag (1 indicates activation of debugging printouts) */
 79:   PetscViewer   viewer1, viewer2;  /* viewers for the solution and error */
 80:   double   norm_2, norm_max;  /* error norms */
 81: } AppCtx;

 83: /* 
 84:    User-defined routines
 85: */
 86: extern int InitialConditions(Vec,AppCtx*);
 87: extern int RHSMatrixHeat(TS,PetscReal,Mat*,Mat*,MatStructure*,void*);
 88: extern int Monitor(TS,int,PetscReal,Vec,void*);
 89: extern int ExactSolution(PetscReal,Vec,AppCtx*);
 90: extern int MyBCRoutine(TS,PetscReal,Vec,void*);

 94: int main(int argc,char **argv)
 95: {
 96:   AppCtx        appctx;                 /* user-defined application context */
 97:   TS            ts;                     /* timestepping context */
 98:   Mat           A;                      /* matrix data structure */
 99:   Vec           u;                      /* approximate solution vector */
100:   double        time_total_max = 100.0; /* default max total time */
101:   int           time_steps_max = 100;   /* default max timesteps */
102:   PetscDraw          draw;                   /* drawing context */
103:   PetscTruth    flg;
104:   int           ierr,  steps, size, m;
105:   double        dt;
106:   PetscReal     ftime;

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:      Initialize program and set problem parameters
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: 
112:   PetscInitialize(&argc,&argv,(char*)0,help);
113:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
114:   if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");

116:   m    = 60;
117:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flg);
118:   PetscOptionsHasName(PETSC_NULL,"-debug",&flg);
119:   appctx.m        = m;
120:   appctx.h        = 1.0/(m-1.0);
121:   appctx.norm_2   = 0.0;
122:   appctx.norm_max = 0.0;
123:   appctx.debug    = (int) flg;
124:   PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");

126:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127:      Create vector data structures
128:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

130:   /* 
131:      Create vector data structures for approximate and exact solutions
132:   */
133:   VecCreateSeq(PETSC_COMM_SELF,m,&u);
134:   VecDuplicate(u,&appctx.solution);

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:      Set up displays to show graphs of the solution and error 
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

140:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
141:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
142:   PetscDrawSetDoubleBuffer(draw);
143:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
144:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
145:   PetscDrawSetDoubleBuffer(draw);

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:      Create timestepping solver context
149:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

151:   TSCreate(PETSC_COMM_SELF,&ts);
152:   TSSetProblemType(ts,TS_LINEAR);

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:      Set optional user-defined monitoring routine
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

162:      Create matrix data structure; set matrix evaluation routine.
163:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

165:   MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,m,m,&A);

167:   PetscOptionsHasName(PETSC_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:     TSSetRHSMatrix(ts,A,A,RHSMatrixHeat,&appctx);
175:   } else {
176:     /*
177:        For linear problems with a time-independent f(u) in the equation 
178:        u_t = f(u), the user provides the discretized right-hand-side
179:        as a matrix only once, and then sets a null matrix evaluation
180:        routine.
181:     */
182:     MatStructure A_structure;
183:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);
184:     TSSetRHSMatrix(ts,A,A,PETSC_NULL,&appctx);
185:   }

187:   /* Treat the problem as having time-dependent boundary conditions */
188:   PetscOptionsHasName(PETSC_NULL,"-time_dependent_bc",&flg);
189:   if (flg) {
190:     TSSetRHSBoundaryConditions(ts,MyBCRoutine,&appctx);
191:   }

193:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194:      Set solution vector and initial timestep
195:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

197:   dt = appctx.h*appctx.h/2.0;
198:   TSSetInitialTimeStep(ts,0.0,dt);
199:   TSSetSolution(ts,u);

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

211:   TSSetDuration(ts,time_steps_max,time_total_max);
212:   TSSetFromOptions(ts);

214:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215:      Solve the problem
216:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

218:   /*
219:      Evaluate initial conditions
220:   */
221:   InitialConditions(u,&appctx);

223:   /*
224:      Run the timestepping solver
225:   */
226:   TSStep(ts,&steps,&ftime);

228:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229:      View timestepping solver info
230:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

232:   PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %g, avg. error (max norm) = %g\n",
233:               appctx.norm_2/steps,appctx.norm_max/steps);
234:   TSView(ts,PETSC_VIEWER_STDOUT_SELF);

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

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

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

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

267:    Output Parameter:
268:    u - vector with solution at initial time (global)
269: */
270: int InitialConditions(Vec u,AppCtx *appctx)
271: {
272:   PetscScalar *u_localptr;
273:   /*PetscScalar h = appctx->h;*/
274:   int    i, ierr;

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

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

297:   /* 
298:      Restore vector
299:   */
300:   VecRestoreArray(u,&u_localptr);

302:   /* 
303:      Print debugging information if desired
304:   */
305:   if (appctx->debug) {
306:      printf("initial guess vector\n");
307:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
308:   }

310:   return 0;
311: }
312: /* --------------------------------------------------------------------- */
315: /*
316:    ExactSolution - Computes the exact solution at a given time.

318:    Input Parameters:
319:    t - current time
320:    solution - vector in which exact solution will be computed
321:    appctx - user-defined application context

323:    Output Parameter:
324:    solution - vector with the newly computed exact solution
325: */
326: int ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
327: {
328:   PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2;
329:   int    i, ierr;

331:   /*
332:      Get a pointer to vector data.
333:   */
334:   VecGetArray(solution,&s_localptr);

336:   /* 
337:      Simply write the solution directly into the array locations.
338:      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
339:   */
340:   ex1 = exp(-36.*PETSC_PI*PETSC_PI*t); ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
341:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
342:   for (i=0; i<appctx->m; i++) {
343:     s_localptr[i] = sin(PetscRealPart(sc1)*(double)i)*ex1 + 3.*sin(PetscRealPart(sc2)*(double)i)*ex2;
344:   }

346:   /* 
347:      Restore vector
348:   */
349:   VecRestoreArray(solution,&s_localptr);
350:   return 0;
351: }
352: /* --------------------------------------------------------------------- */
355: /*
356:    Monitor - User-provided routine to monitor the solution computed at 
357:    each timestep.  This example plots the solution and computes the
358:    error in two different norms.

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

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

381:   /* 
382:      View a graph of the current iterate
383:   */
384:   VecView(u,appctx->viewer2);

386:   /* 
387:      Compute the exact solution
388:   */
389:   ExactSolution(ctime,appctx->solution,appctx);

391:   /*
392:      Print debugging information if desired
393:   */
394:   if (appctx->debug) {
395:      printf("Computed solution vector\n");
396:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
397:      printf("Exact solution vector\n");
398:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
399:   }

401:   /*
402:      Compute the 2-norm and max-norm of the error
403:   */
404:   VecAXPY(&mone,u,appctx->solution);
405:   VecNorm(appctx->solution,NORM_2,&norm_2);
406:   norm_2 = sqrt(appctx->h)*norm_2;
407:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

409:   TSGetTimeStep(ts,&dt);
410:   printf("Timestep %d: step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n",
411:          step,dt,ctime,norm_2,norm_max);
412:   appctx->norm_2   += norm_2;
413:   appctx->norm_max += norm_max;

415:   dttol = .0001;
416:   PetscOptionsGetReal(PETSC_NULL,"-dttol",&dttol,&flg);
417:   if (dt < dttol) {
418:     dt *= .999;
419:     TSSetTimeStep(ts,dt);
420:   }

422:   /* 
423:      View a graph of the error
424:   */
425:   VecView(appctx->solution,appctx->viewer1);

427:   /*
428:      Print debugging information if desired
429:   */
430:   if (appctx->debug) {
431:      printf("Error vector\n");
432:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
433:   }

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

444:    Input Parameters:
445:    ts - the TS context
446:    t - current time
447:    global_in - global input vector
448:    dummy - optional user-defined context, as set by TSetRHSJacobian()

450:    Output Parameters:
451:    AA - Jacobian matrix
452:    BB - optionally different preconditioning matrix
453:    str - flag indicating matrix structure

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

468:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
469:      Compute entries for the locally owned part of the matrix
470:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
471:   /* 
472:      Set matrix rows corresponding to boundary data
473:   */

475:   mstart = 0;
476:   v[0] = 1.0;
477:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
478:   mstart++;

480:   mend--;
481:   v[0] = 1.0;
482:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

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

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

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

525:   /*
526:      Set and option to indicate that we will never add a new nonzero location 
527:      to the matrix. If we do, it will generate an error.
528:   */
529:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);

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

549:   VecGetArray(f,&fa);
550:   fa[0] = 0.0;
551:   fa[m-1] = 1.0;
552:   VecRestoreArray(f,&fa);
553:   printf("t=%g\n",t);
554: 
555:   return 0;
556: }