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: #undef __FUNCT__
 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:   PetscTruth    flg;
 97:   PetscReal     dt,ftime;

 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:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179:      Set solution vector and initial timestep
180:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

196:   TSSetDuration(ts,time_steps_max,time_total_max);
197:   TSSetFromOptions(ts);

199:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200:      Solve the problem
201:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

213:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214:      View timestepping solver info
215:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

248:    Input Parameter:
249:    u - uninitialized solution vector (global)
250:    appctx - user-defined application context

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

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

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

280:   /* 
281:      Restore vector
282:   */
283:   VecRestoreArray(u,&u_localptr);

285:   /* 
286:      Print debugging information if desired
287:   */
288:   if (appctx->debug) {
289:      printf("initial guess vectorn");
290:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
291:   }

293:   return 0;
294: }
295: /* --------------------------------------------------------------------- */
296: #undef __FUNCT__
298: /*
299:    ExactSolution - Computes the exact solution at a given time.

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

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

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

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

329:   /* 
330:      Restore vector
331:   */
332:   VecRestoreArray(solution,&s_localptr);
333:   return 0;
334: }
335: /* --------------------------------------------------------------------- */
336: #undef __FUNCT__
338: /*
339:    Monitor - User-provided routine to monitor the solution computed at 
340:    each timestep.  This example plots the solution and computes the
341:    error in two different norms.

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

361:   /* 
362:      View a graph of the current iterate
363:   */
364:   VecView(u,appctx->viewer2);

366:   /* 
367:      Compute the exact solution
368:   */
369:   ExactSolution(time,appctx->solution,appctx);

371:   /*
372:      Print debugging information if desired
373:   */
374:   if (appctx->debug) {
375:      printf("Computed solution vectorn");
376:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
377:      printf("Exact solution vectorn");
378:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
379:   }

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

389:   printf("Timestep %d: time = %g, 2-norm error = %g, max norm error = %gn",
390:          step,time,norm_2,norm_max);
391:   appctx->norm_2   += norm_2;
392:   appctx->norm_max += norm_max;

394:   /* 
395:      View a graph of the error
396:   */
397:   VecView(appctx->solution,appctx->viewer1);

399:   /*
400:      Print debugging information if desired
401:   */
402:   if (appctx->debug) {
403:      printf("Error vectorn");
404:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
405:   }

407:   return 0;
408: }
409: /* --------------------------------------------------------------------- */
410: #undef __FUNCT__
412: /*
413:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
414:    matrix for the heat equation.

416:    Input Parameters:
417:    ts - the TS context
418:    t - current time
419:    global_in - global input vector
420:    dummy - optional user-defined context, as set by TSetRHSJacobian()

422:    Output Parameters:
423:    AA - Jacobian matrix
424:    BB - optionally different preconditioning matrix
425:    str - flag indicating matrix structure

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

440:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
441:      Compute entries for the locally owned part of the matrix
442:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
443:   /* 
444:      Set matrix rows corresponding to boundary data
445:   */

447:   mstart = 0;
448:   v[0] = 1.0;
449:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
450:   mstart++;

452:   mend--;
453:   v[0] = 1.0;
454:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

456:   /*
457:      Set matrix rows corresponding to interior data.  We construct the 
458:      matrix one row at a time.
459:   */
460:   v[0] = sone; v[1] = stwo; v[2] = sone;
461:   for (i=mstart; i<mend; i++) {
462:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
463:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
464:   }

466:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
467:      Complete the matrix assembly process and set some options
468:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
469:   /*
470:      Assemble matrix, using the 2-step process:
471:        MatAssemblyBegin(), MatAssemblyEnd()
472:      Computations can be done while messages are in transition
473:      by placing code between these two statements.
474:   */
475:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
476:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

497:   /*
498:      Set and option to indicate that we will never add a new nonzero location 
499:      to the matrix. If we do, it will generate an error.
500:   */
501:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);

503:   return 0;
504: }