Actual source code: ex3.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:    Processors: 1
 16: */

 18: /* ------------------------------------------------------------------------

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

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

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

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

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

 48:   ------------------------------------------------------------------------- */

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

 60: #include <petscts.h>

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

 75: /* 
 76:    User-defined routines
 77: */
 78: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
 79: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
 80: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 81: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
 82: extern PetscErrorCode 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:   PetscInt       time_steps_max = 100;   /* default max timesteps */
 94:   PetscDraw      draw;                   /* drawing context */
 96:   PetscInt       steps,m;
 97:   PetscMPIInt    size;
 98:   PetscReal      dt;
 99:   PetscBool      flg;

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

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

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

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

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

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

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

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

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

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

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

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

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

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

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

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_final_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:   TSSolve(ts,u,PETSC_NULL);
218:   TSGetTimeStepNumber(ts,&steps);

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

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

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

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

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

255:    Input Parameter:
256:    u - uninitialized solution vector (global)
257:    appctx - user-defined application context

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

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

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

288:   /* 
289:      Restore vector
290:   */
291:   VecRestoreArray(u,&u_localptr);

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

301:   return 0;
302: }
303: /* --------------------------------------------------------------------- */
306: /*
307:    ExactSolution - Computes the exact solution at a given time.

309:    Input Parameters:
310:    t - current time
311:    solution - vector in which exact solution will be computed
312:    appctx - user-defined application context

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

323:   /*
324:      Get a pointer to vector data.
325:   */
326:   VecGetArray(solution,&s_localptr);

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

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

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

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

376:   /* 
377:      Compute the exact solution
378:   */
379:   ExactSolution(time,appctx->solution,appctx);

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

391:   /*
392:      Compute the 2-norm and max-norm of the error
393:   */
394:   VecAXPY(appctx->solution,-1.0,u);
395:   VecNorm(appctx->solution,NORM_2,&norm_2);
396:   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
397:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

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

405:   dttol = .0001;
406:   PetscOptionsGetReal(PETSC_NULL,"-dttol",&dttol,PETSC_NULL);
407:   if (dt < dttol) {
408:     dt *= .999;
409:     TSSetTimeStep(ts,dt);
410:   }

412:   /* 
413:      View a graph of the error
414:   */
415:   VecView(appctx->solution,appctx->viewer1);

417:   /*
418:      Print debugging information if desired
419:   */
420:   if (appctx->debug) {
421:      printf("Error vector\n");
422:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
423:   }

425:   return 0;
426: }
427: /* --------------------------------------------------------------------- */
430: /*
431:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
432:    matrix for the heat equation.

434:    Input Parameters:
435:    ts - the TS context
436:    t - current time
437:    global_in - global input vector
438:    dummy - optional user-defined context, as set by TSetRHSJacobian()

440:    Output Parameters:
441:    AA - Jacobian matrix
442:    BB - optionally different preconditioning matrix
443:    str - flag indicating matrix structure

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

459:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
460:      Compute entries for the locally owned part of the matrix
461:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
462:   /* 
463:      Set matrix rows corresponding to boundary data
464:   */

466:   mstart = 0;
467:   v[0] = 1.0;
468:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
469:   mstart++;

471:   mend--;
472:   v[0] = 1.0;
473:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

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

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

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

516:   /*
517:      Set and option to indicate that we will never add a new nonzero location 
518:      to the matrix. If we do, it will generate an error.
519:   */
520:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

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

540:   VecGetArray(f,&fa);
541:   fa[0] = 0.0;
542:   fa[m-1] = 0.0;
543:   VecRestoreArray(f,&fa);
544:   printf("t=%g\n",t);
545: 
546:   return 0;
547: }