Actual source code: ex6.c

petsc-3.3-p7 2013-05-11
  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:   -debug              : Activate debugging printouts\n\
  9:   -nox                : Deactivate x-window graphics\n\n";

 11: /*
 12:    Concepts: TS^time-dependent linear problems
 13:    Concepts: TS^heat equation
 14:    Concepts: TS^diffusion equation
 15:    Routines: TSCreate(); TSSetSolution(); TSSetRHSJacobian(), TSSetIJacobian();
 16:    Routines: TSSetInitialTimeStep(); TSSetDuration(); TSMonitorSet();
 17:    Routines: TSSetFromOptions(); TSStep(); TSDestroy(); 
 18:    Routines: TSSetTimeStep(); TSGetTimeStep();
 19:    Processors: 1
 20: */

 22: /* ------------------------------------------------------------------------

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

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

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

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

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

 52:   ------------------------------------------------------------------------- */

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

 64: #include <petscts.h>

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

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

 90: int main(int argc,char **argv)
 91: {
 92:   AppCtx         appctx;                 /* user-defined application context */
 93:   TS             ts;                     /* timestepping context */
 94:   Mat            A;                      /* matrix data structure */
 95:   Vec            u;                      /* approximate solution vector */
 96:   PetscReal      time_total_max = 100.0; /* default max total time */
 97:   PetscInt       time_steps_max = 100;   /* default max timesteps */
 98:   PetscDraw      draw;                   /* drawing context */
100:   PetscInt       steps, m;
101:   PetscMPIInt    size;
102:   PetscReal      dt;
103:   PetscReal      ftime;
104:   PetscBool      flg;
105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:      Initialize program and set problem parameters
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: 
109:   PetscInitialize(&argc,&argv,(char*)0,help);
110:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
111:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

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

122:   PetscOptionsGetInt(PETSC_NULL,"-time_steps_max",&time_steps_max,PETSC_NULL);

124:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125:      Create vector data structures
126:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:      Create timestepping solver context
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

149:   TSCreate(PETSC_COMM_SELF,&ts);
150:   TSSetProblemType(ts,TS_LINEAR);

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:      Set optional user-defined monitoring routine
154:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

156:   TSMonitorSet(ts,Monitor,&appctx,PETSC_NULL);

158:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

160:      Create matrix data structure; set matrix evaluation routine.
161:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

163:   MatCreate(PETSC_COMM_SELF,&A);
164:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
165:   MatSetFromOptions(A);
166:   MatSetUp(A);

168:   PetscOptionsHasName(PETSC_NULL,"-time_dependent_rhs",&flg);
169:   if (flg) {
170:     /*
171:        For linear problems with a time-dependent f(u,t) in the equation 
172:        u_t = f(u,t), the user provides the discretized right-hand-side
173:        as a time-dependent matrix.
174:     */
175:     TSSetRHSFunction(ts,PETSC_NULL,TSComputeRHSFunctionLinear,&appctx);
176:     TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
177:   } else {
178:     /*
179:        For linear problems with a time-independent f(u) in the equation 
180:        u_t = f(u), the user provides the discretized right-hand-side
181:        as a matrix only once, and then sets a null matrix evaluation
182:        routine.
183:     */
184:     MatStructure A_structure;
185:     RHSMatrixHeat(ts,0.0,u,&A,&A,&A_structure,&appctx);
186:     TSSetRHSFunction(ts,PETSC_NULL,TSComputeRHSFunctionLinear,&appctx);
187:     TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
188:   }

190:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191:      Set solution vector and initial timestep
192:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

194:   dt = appctx.h*appctx.h/2.0;
195:   TSSetInitialTimeStep(ts,0.0,dt);
196:   TSSetSolution(ts,u);

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

208:   TSSetDuration(ts,time_steps_max,time_total_max);
209:   TSSetFromOptions(ts);

211:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212:      Solve the problem
213:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

220:   /*
221:      Run the timestepping solver
222:   */
223:   TSSolve(ts,u,&ftime);
224:   TSGetTimeStepNumber(ts,&steps);

226:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227:      View timestepping solver info
228:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

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

265:    Output Parameter:
266:    u - vector with solution at initial time (global)
267: */
268: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
269: {
270:   PetscScalar    *u_localptr;
271:   PetscInt       i;

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

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

294:   /* 
295:      Restore vector
296:   */
297:   VecRestoreArray(u,&u_localptr);

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

307:   return 0;
308: }
309: /* --------------------------------------------------------------------- */
312: /*
313:    ExactSolution - Computes the exact solution at a given time.

315:    Input Parameters:
316:    t - current time
317:    solution - vector in which exact solution will be computed
318:    appctx - user-defined application context

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

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

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

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

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

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

378:   /* 
379:      View a graph of the current iterate
380:   */
381:   VecView(u,appctx->viewer2);

383:   /* 
384:      Compute the exact solution
385:   */
386:   ExactSolution(crtime,appctx->solution,appctx);

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

398:   /*
399:      Compute the 2-norm and max-norm of the error
400:   */
401:   VecAXPY(appctx->solution,-1.0,u);
402:   VecNorm(appctx->solution,NORM_2,&norm_2);
403:   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
404:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

406:   TSGetTimeStep(ts,&dt);
407:   if (norm_2 > 1.e-2){
408:     printf("Timestep %d: step size = %G, time = %G, 2-norm error = %G, max norm error = %G\n",
409:          (int)step,dt,crtime,norm_2,norm_max);
410:   }
411:   appctx->norm_2   += norm_2;
412:   appctx->norm_max += norm_max;

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

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

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

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

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

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

454:    Notes:
455:    Recall that MatSetValues() uses 0-based row and column numbers
456:    in Fortran as well as in C.
457: */
458: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
459: {
460:   Mat            A = *AA;                      /* Jacobian matrix */
461:   AppCtx         *appctx = (AppCtx *) ctx;     /* user-defined application context */
462:   PetscInt       mstart = 0;
463:   PetscInt       mend = appctx->m;
465:   PetscInt       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,PETSC_TRUE);

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: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
544: {
545:   AppCtx         *appctx = (AppCtx *) ctx;     /* user-defined application context */
547:   PetscInt       m = appctx->m;
548:   PetscScalar    *fa;

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