Actual source code: ex5.c

  1: /*$Id: ex5.c,v 1.25 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:   -debug              : Activate debugging printoutsn
 10:   -nox                : Deactivate x-window graphicsnn";

 12: /*
 13:    Concepts: TS^time-dependent linear problems
 14:    Concepts: TS^heat equation
 15:    Concepts: TS^diffusion equation
 16:    Processors: 1
 17: */

 19: /* ------------------------------------------------------------------------

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

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

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

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

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

 49:   ------------------------------------------------------------------------- */

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

 61:  #include petscts.h

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

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

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

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

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

114:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115:      Create vector data structures
116:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

118:   /* 
119:      Create vector data structures for approximate and exact solutions
120:   */
121:   VecCreateSeq(PETSC_COMM_SELF,m,&u);
122:   VecDuplicate(u,&appctx.solution);

124:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125:      Set up displays to show graphs of the solution and error 
126:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

128:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);
129:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
130:   PetscDrawSetDoubleBuffer(draw);
131:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);
132:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
133:   PetscDrawSetDoubleBuffer(draw);

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:      Create timestepping solver context
137:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

139:   TSCreate(PETSC_COMM_SELF,&ts);
140:   TSSetProblemType(ts,TS_LINEAR);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:      Set optional user-defined monitoring routine
144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

150:      Create matrix data structure; set matrix evaluation routine.
151:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

153:   MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,m,m,&A);
154:   MatSetFromOptions(A);

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

176:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177:      Set solution vector and initial timestep
178:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

180:   dt = appctx.h*appctx.h/2.0;
181:   TSSetInitialTimeStep(ts,0.0,dt);
182:   TSSetSolution(ts,u);

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:      Customize timestepping solver:  
186:        - Set the solution method to be the Backward Euler method.
187:        - Set timestepping duration info 
188:      Then set runtime options, which can override these defaults.
189:      For example,
190:           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
191:      to override the defaults set by TSSetDuration().
192:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

194:   TSSetDuration(ts,time_steps_max,time_total_max);
195:   TSSetFromOptions(ts);

197:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198:      Solve the problem
199:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

201:   /*
202:      Evaluate initial conditions
203:   */
204:   InitialConditions(u,&appctx);

206:   /*
207:      Run the timestepping solver
208:   */
209:   TSStep(ts,&steps,&ftime);

211:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212:      View timestepping solver info
213:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

219:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220:      Free work space.  All PETSc objects should be destroyed when they
221:      are no longer needed.
222:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

224:   TSDestroy(ts);
225:   MatDestroy(A);
226:   VecDestroy(u);
227:   PetscViewerDestroy(appctx.viewer1);
228:   PetscViewerDestroy(appctx.viewer2);
229:   VecDestroy(appctx.solution);

231:   /*
232:      Always call PetscFinalize() before exiting a program.  This routine
233:        - finalizes the PETSc libraries as well as MPI
234:        - provides summary and diagnostic information if certain runtime
235:          options are chosen (e.g., -log_summary). 
236:   */
237:   PetscFinalize();
238:   return 0;
239: }
240: /* --------------------------------------------------------------------- */
241: /*
242:    InitialConditions - Computes the solution at the initial time. 

244:    Input Parameter:
245:    u - uninitialized solution vector (global)
246:    appctx - user-defined application context

248:    Output Parameter:
249:    u - vector with solution at initial time (global)
250: */
251: int InitialConditions(Vec u,AppCtx *appctx)
252: {
253:   PetscScalar *u_localptr,h = appctx->h;
254:   int    i,ierr;

256:   /* 
257:     Get a pointer to vector data.
258:     - For default PETSc vectors, VecGetArray() returns a pointer to
259:       the data array.  Otherwise, the routine is implementation dependent.
260:     - You MUST call VecRestoreArray() when you no longer need access to
261:       the array.
262:     - Note that the Fortran interface to VecGetArray() differs from the
263:       C version.  See the users manual for details.
264:   */
265:   VecGetArray(u,&u_localptr);

267:   /* 
268:      We initialize the solution array by simply writing the solution
269:      directly into the array locations.  Alternatively, we could use
270:      VecSetValues() or VecSetValuesLocal().
271:   */
272:   for (i=0; i<appctx->m; i++) {
273:     u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h);
274:   }

276:   /* 
277:      Restore vector
278:   */
279:   VecRestoreArray(u,&u_localptr);

281:   /* 
282:      Print debugging information if desired
283:   */
284:   if (appctx->debug) {
285:      printf("initial guess vectorn");
286:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
287:   }

289:   return 0;
290: }
291: /* --------------------------------------------------------------------- */
292: /*
293:    ExactSolution - Computes the exact solution at a given time.

295:    Input Parameters:
296:    t - current time
297:    solution - vector in which exact solution will be computed
298:    appctx - user-defined application context

300:    Output Parameter:
301:    solution - vector with the newly computed exact solution
302: */
303: int ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
304: {
305:   PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
306:   int         i,ierr;

308:   /*
309:      Get a pointer to vector data.
310:   */
311:   VecGetArray(solution,&s_localptr);

313:   /* 
314:      Simply write the solution directly into the array locations.
315:      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
316:   */
317:   ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
318:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
319:   for (i=0; i<appctx->m; i++) {
320:     s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2;
321:   }

323:   /* 
324:      Restore vector
325:   */
326:   VecRestoreArray(solution,&s_localptr);
327:   return 0;
328: }
329: /* --------------------------------------------------------------------- */
330: /*
331:    Monitor - User-provided routine to monitor the solution computed at 
332:    each timestep.  This example plots the solution and computes the
333:    error in two different norms.

335:    Input Parameters:
336:    ts     - the timestep context
337:    step   - the count of the current step (with 0 meaning the
338:              initial condition)
339:    time   - the current time
340:    u      - the solution at this timestep
341:    ctx    - the user-provided context for this monitoring routine.
342:             In this case we use the application context which contains 
343:             information about the problem size, workspace and the exact 
344:             solution.
345: */
346: int Monitor(TS ts,int step,PetscReal time,Vec u,void *ctx)
347: {
348:   AppCtx      *appctx = (AppCtx*) ctx;   /* user-defined application context */
349:   int         ierr;
350:   PetscReal   norm_2,norm_max;
351:   PetscScalar mone = -1.0;

353:   /* 
354:      View a graph of the current iterate
355:   */
356:   VecView(u,appctx->viewer2);

358:   /* 
359:      Compute the exact solution
360:   */
361:   ExactSolution(time,appctx->solution,appctx);

363:   /*
364:      Print debugging information if desired
365:   */
366:   if (appctx->debug) {
367:      printf("Computed solution vectorn");
368:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
369:      printf("Exact solution vectorn");
370:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
371:   }

373:   /*
374:      Compute the 2-norm and max-norm of the error
375:   */
376:   VecAXPY(&mone,u,appctx->solution);
377:   VecNorm(appctx->solution,NORM_2,&norm_2);
378:   norm_2 = sqrt(appctx->h)*norm_2;
379:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

381:   printf("Timestep %d: time = %g, 2-norm error = %g, max norm error = %gn",
382:          step,time,norm_2,norm_max);
383:   appctx->norm_2   += norm_2;
384:   appctx->norm_max += norm_max;

386:   /* 
387:      View a graph of the error
388:   */
389:   VecView(appctx->solution,appctx->viewer1);

391:   /*
392:      Print debugging information if desired
393:   */
394:   if (appctx->debug) {
395:      printf("Error vectorn");
396:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
397:   }

399:   return 0;
400: }
401: /* --------------------------------------------------------------------- */
402: /*
403:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
404:    matrix for the heat equation.

406:    Input Parameters:
407:    ts - the TS context
408:    t - current time
409:    global_in - global input vector
410:    dummy - optional user-defined context, as set by TSetRHSJacobian()

412:    Output Parameters:
413:    AA - Jacobian matrix
414:    BB - optionally different preconditioning matrix
415:    str - flag indicating matrix structure

417:   Notes:
418:   Recall that MatSetValues() uses 0-based row and column numbers
419:   in Fortran as well as in C.
420: */
421: int RHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
422: {
423:   Mat    A = *AA;                      /* Jacobian matrix */
424:   AppCtx *appctx = (AppCtx*)ctx;     /* user-defined application context */
425:   int    mstart = 0;
426:   int    mend = appctx->m;
427:   int    ierr,i,idx[3];
428:   PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;

430:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
431:      Compute entries for the locally owned part of the matrix
432:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
433:   /* 
434:      Set matrix rows corresponding to boundary data
435:   */

437:   mstart = 0;
438:   v[0] = 1.0;
439:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
440:   mstart++;

442:   mend--;
443:   v[0] = 1.0;
444:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

446:   /*
447:      Set matrix rows corresponding to interior data.  We construct the 
448:      matrix one row at a time.
449:   */
450:   v[0] = sone; v[1] = stwo; v[2] = sone;
451:   for (i=mstart; i<mend; i++) {
452:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
453:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
454:   }

456:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
457:      Complete the matrix assembly process and set some options
458:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
459:   /*
460:      Assemble matrix, using the 2-step process:
461:        MatAssemblyBegin(), MatAssemblyEnd()
462:      Computations can be done while messages are in transition
463:      by placing code between these two statements.
464:   */
465:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
466:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

468:   /*
469:      Set flag to indicate that the Jacobian matrix retains an identical
470:      nonzero structure throughout all timestepping iterations (although the
471:      values of the entries change). Thus, we can save some work in setting
472:      up the preconditioner (e.g., no need to redo symbolic factorization for
473:      ILU/ICC preconditioners).
474:       - If the nonzero structure of the matrix is different during
475:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
476:         must be used instead.  If you are unsure whether the matrix
477:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
478:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
479:         believes your assertion and does not check the structure
480:         of the matrix.  If you erroneously claim that the structure
481:         is the same when it actually is not, the new preconditioner
482:         will not function correctly.  Thus, use this optimization
483:         feature with caution!
484:   */
485:   *str = SAME_NONZERO_PATTERN;

487:   /*
488:      Set and option to indicate that we will never add a new nonzero location 
489:      to the matrix. If we do, it will generate an error.
490:   */
491:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);

493:   return 0;
494: }