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: #undef __FUNCT__
 88: int main(int argc,char **argv)
 89: {
 90:   AppCtx        appctx;                 /* user-defined application context */
 91:   TS            ts;                     /* timestepping context */
 92:   Mat           A;                      /* matrix data structure */
 93:   Vec           u;                      /* approximate solution vector */
 94:   PetscReal     time_total_max = 100.0; /* default max total time */
 95:   int           time_steps_max = 100;   /* default max timesteps */
 96:   PetscDraw     draw;                   /* drawing context */
 97:   int           ierr,steps,size,m;
 98:   PetscReal     dt,ftime;
 99:   PetscTruth    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(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 processorn");

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

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

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

157:   MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,m,m,&A);
158:   MatSetFromOptions(A);

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

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

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

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

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

204:   TSSetDuration(ts,time_steps_max,time_total_max);
205:   TSSetFromOptions(ts);

207:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208:      Solve the problem
209:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

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

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

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

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

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

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

301:   return 0;
302: }
303: /* --------------------------------------------------------------------- */
304: #undef __FUNCT__
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: int ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
318: {
319:   PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
320:   int         i,ierr;

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

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

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

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

354:    Input Parameters:
355:    ts     - the timestep context
356:    step   - the count of the current step (with 0 meaning the
357:              initial condition)
358:    time   - the current time
359:    u      - the solution at this timestep
360:    ctx    - the user-provided context for this monitoring routine.
361:             In this case we use the application context which contains 
362:             information about the problem size, workspace and the exact 
363:             solution.
364: */
365: int Monitor(TS ts,int step,PetscReal time,Vec u,void *ctx)
366: {
367:   AppCtx       *appctx = (AppCtx*) ctx;   /* user-defined application context */
368:   int          ierr;
369:   PetscReal    norm_2,norm_max,dt,dttol;
370:   PetscScalar  mone = -1.0;
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 vectorn");
386:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
387:      printf("Exact solution vectorn");
388:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
389:   }

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

399:   TSGetTimeStep(ts,&dt);
400:   printf("Timestep %3d: step size = %-11g, time = %-11g, 2-norm error = %-11g, max norm error = %-11gn",
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 vectorn");
422:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
423:   }

425:   return 0;
426: }
427: /* --------------------------------------------------------------------- */
428: #undef __FUNCT__
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: int RHSMatrixHeat(TS ts,PetscReal t,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:   int    mstart = 0;
454:   int    mend = appctx->m;
455:   int    ierr,i,idx[3];
456:   PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;

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

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

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

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

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

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

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

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

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