Actual source code: ex3.c

  1: /*$Id: ex3.c,v 1.28 2001/08/10 03:34:17 bsmith Exp $*/

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

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

 13: /*
 14:    Concepts: TS^time-dependent linear problems
 15:    Concepts: TS^heat equation
 16:    Concepts: TS^diffusion equation
 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 "petscts.h" so that we can use TS solvers.  Note that this file
 54:    automatically includes:
 55:      petsc.h       - base PETSc routines   petscvec.h  - vectors
 56:      petscsys.h    - system routines       petscmat.h  - matrices
 57:      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
 58:      petscviewer.h - viewers               petscpc.h   - preconditioners
 59:      petscsles.h   - linear solvers        petscsnes.h - nonlinear solvers
 60: */

 62:  #include petscts.h

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

 77: /* 
 78:    User-defined routines
 79: */
 80: extern int InitialConditions(Vec,AppCtx*);
 81: extern int RHSMatrixHeat(TS,PetscReal,Mat*,Mat*,MatStructure*,void*);
 82: extern int Monitor(TS,int,PetscReal,Vec,void*);
 83: extern int ExactSolution(PetscReal,Vec,AppCtx*);
 84: extern int MyBCRoutine(TS,PetscReal,Vec,void*);

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

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

107:   m    = 60;
108:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
109:   PetscOptionsHasName(PETSC_NULL,"-debug",&appctx.debug);
110:   appctx.m        = m;
111:   appctx.h        = 1.0/(m-1.0);
112:   appctx.norm_2   = 0.0;
113:   appctx.norm_max = 0.0;
114:   PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processorn");

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:   TSSetMonitor(ts,Monitor,&appctx,PETSC_NULL);

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

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

155:   MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,m,m,&A);
156:   MatSetFromOptions(A);

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

178:   /* Treat the problem as having time-dependent boundary conditions */
179:   PetscOptionsHasName(PETSC_NULL,"-time_dependent_bc",&flg);
180:   if (flg) {
181:     TSSetRHSBoundaryConditions(ts,MyBCRoutine,&appctx);
182:   }

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:      Set solution vector and initial timestep
186:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

202:   TSSetDuration(ts,time_steps_max,time_total_max);
203:   TSSetFromOptions(ts);

205:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206:      Solve the problem
207:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

214:   /*
215:      Run the timestepping solver
216:   */
217:   TSStep(ts,&steps,&ftime);

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

223:   PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %g, avg. error (max norm) = %gn",
224:               appctx.norm_2/steps,appctx.norm_max/steps);
225:   TSView(ts,PETSC_VIEWER_STDOUT_SELF);

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

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

239:   /*
240:      Always call PetscFinalize() before exiting a program.  This routine
241:        - finalizes the PETSc libraries as well as MPI
242:        - provides summary and diagnostic information if certain runtime
243:          options are chosen (e.g., -log_summary). 
244:   */
245:   PetscFinalize();
246:   return 0;
247: }
248: /* --------------------------------------------------------------------- */
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: int InitialConditions(Vec u,AppCtx *appctx)
260: {
261:   PetscScalar *u_localptr,h = appctx->h;
262:   int    i,ierr;

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

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

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 vectorn");
294:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
295:   }

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

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

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

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

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

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:    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: int Monitor(TS ts,int step,PetscReal time,Vec u,void *ctx)
358: {
359:   AppCtx       *appctx = (AppCtx*) ctx;   /* user-defined application context */
360:   int          ierr;
361:   PetscReal    norm_2,norm_max,dt,dttol;
362:   PetscScalar  mone = -1.0;
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 vectorn");
378:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
379:      printf("Exact solution vectorn");
380:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
381:   }

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

391:   TSGetTimeStep(ts,&dt);
392:   printf("Timestep %3d: step size = %-11g, time = %-11g, 2-norm error = %-11g, max norm error = %-11gn",
393:          step,dt,time,norm_2,norm_max);
394:   appctx->norm_2   += norm_2;
395:   appctx->norm_max += norm_max;

397:   dttol = .0001;
398:   PetscOptionsGetReal(PETSC_NULL,"-dttol",&dttol,PETSC_NULL);
399:   if (dt < dttol) {
400:     dt *= .999;
401:     TSSetTimeStep(ts,dt);
402:   }

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

409:   /*
410:      Print debugging information if desired
411:   */
412:   if (appctx->debug) {
413:      printf("Error vectorn");
414:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
415:   }

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

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

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

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

448:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
449:      Compute entries for the locally owned part of the matrix
450:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
451:   /* 
452:      Set matrix rows corresponding to boundary data
453:   */

455:   mstart = 0;
456:   v[0] = 1.0;
457:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
458:   mstart++;

460:   mend--;
461:   v[0] = 1.0;
462:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

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

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

486:   /*
487:      Set flag to indicate that the Jacobian matrix retains an identical
488:      nonzero structure throughout all timestepping iterations (although the
489:      values of the entries change). Thus, we can save some work in setting
490:      up the preconditioner (e.g., no need to redo symbolic factorization for
491:      ILU/ICC preconditioners).
492:       - If the nonzero structure of the matrix is different during
493:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
494:         must be used instead.  If you are unsure whether the matrix
495:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
496:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
497:         believes your assertion and does not check the structure
498:         of the matrix.  If you erroneously claim that the structure
499:         is the same when it actually is not, the new preconditioner
500:         will not function correctly.  Thus, use this optimization
501:         feature with caution!
502:   */
503:   *str = SAME_NONZERO_PATTERN;

505:   /*
506:      Set and option to indicate that we will never add a new nonzero location 
507:      to the matrix. If we do, it will generate an error.
508:   */
509:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);

511:   return 0;
512: }
513: /* --------------------------------------------------------------------- */
514: /*
515:    Input Parameters:
516:    ts - the TS context
517:    t - current time
518:    f - function
519:    ctx - optional user-defined context, as set by TSetBCFunction()
520:  */
521: int MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
522: {
523:   AppCtx *appctx = (AppCtx*)ctx;     /* user-defined application context */
524:   int    ierr,m = appctx->m;
525:   PetscScalar *fa;

527:   VecGetArray(f,&fa);
528:   fa[0] = 0.0;
529:   fa[m-1] = 0.0;
530:   VecRestoreArray(f,&fa);
531:   printf("t=%gn",t);
532: 
533:   return 0;
534: }