Actual source code: ex6.c

  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:   -time_dependent_bc  : Treat the problem as having time-dependent boundary conditions\n\
  9:   -debug              : Activate debugging printouts\n\
 10:   -nox                : Deactivate x-window graphics\n\n";

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

 23: /* ------------------------------------------------------------------------

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

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

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

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

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

 53:   ------------------------------------------------------------------------- */

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

 65:  #include petscts.h

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

 80: /* 
 81:    User-defined routines
 82: */

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

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

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

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

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

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

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

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

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

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

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

159:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

164:   MatCreate(PETSC_COMM_SELF,&A);
165:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
166:   MatSetFromOptions(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:     TSSetRHSMatrix(ts,A,A,RHSMatrixHeat,&appctx);
176:   } else {
177:     /*
178:        For linear problems with a time-independent f(u) in the equation 
179:        u_t = f(u), the user provides the discretized right-hand-side
180:        as a matrix only once, and then sets a null matrix evaluation
181:        routine.
182:     */
183:     MatStructure A_structure;
184:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);
185:     TSSetRHSMatrix(ts,A,A,PETSC_NULL,&appctx);
186:   }

188:   /* Treat the problem as having time-dependent boundary conditions */
189:   PetscOptionsHasName(PETSC_NULL,"-time_dependent_bc",&flg);
190:   if (flg) {
191:     TSSetRHSBoundaryConditions(ts,MyBCRoutine,&appctx);
192:   }

194:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195:      Set solution vector and initial timestep
196:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

198:   dt = appctx.h*appctx.h/2.0;
199:   TSSetInitialTimeStep(ts,0.0,dt);
200:   TSSetSolution(ts,u);

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

212:   TSSetDuration(ts,time_steps_max,time_total_max);
213:   TSSetFromOptions(ts);

215:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216:      Solve the problem
217:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

219:   /*
220:      Evaluate initial conditions
221:   */
222:   InitialConditions(u,&appctx);

224:   /*
225:      Run the timestepping solver
226:   */
227:   TSStep(ts,&steps,&ftime);

229:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230:      View timestepping solver info
231:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

237:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238:      Free work space.  All PETSc objects should be destroyed when they
239:      are no longer needed.
240:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

242:   TSDestroy(ts);
243:   MatDestroy(A);
244:   VecDestroy(u);
245:   PetscViewerDestroy(appctx.viewer1);
246:   PetscViewerDestroy(appctx.viewer2);
247:   VecDestroy(appctx.solution);

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

264:    Input Parameter:
265:    u - uninitialized solution vector (global)
266:    appctx - user-defined application context

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

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

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

297:   /* 
298:      Restore vector
299:   */
300:   VecRestoreArray(u,&u_localptr);

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

310:   return 0;
311: }
312: /* --------------------------------------------------------------------- */
315: /*
316:    ExactSolution - Computes the exact solution at a given time.

318:    Input Parameters:
319:    t - current time
320:    solution - vector in which exact solution will be computed
321:    appctx - user-defined application context

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

332:   /*
333:      Get a pointer to vector data.
334:   */
335:   VecGetArray(solution,&s_localptr);

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

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

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

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

381:   /* 
382:      View a graph of the current iterate
383:   */
384:   VecView(u,appctx->viewer2);

386:   /* 
387:      Compute the exact solution
388:   */
389:   ExactSolution(crtime,appctx->solution,appctx);

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

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

409:   TSGetTimeStep(ts,&dt);
410:   printf("Timestep %d: step size = %G, time = %G, 2-norm error = %G, max norm error = %G\n",
411:          (int)step,dt,crtime,norm_2,norm_max);
412:   appctx->norm_2   += norm_2;
413:   appctx->norm_max += norm_max;

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

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

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

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

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

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

455:    Notes:
456:    Recall that MatSetValues() uses 0-based row and column numbers
457:    in Fortran as well as in C.
458: */
459: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
460: {
461:   Mat            A = *AA;                      /* Jacobian matrix */
462:   AppCtx         *appctx = (AppCtx *) ctx;     /* user-defined application context */
463:   PetscInt       mstart = 0;
464:   PetscInt       mend = appctx->m;
466:   PetscInt       i, idx[3];
467:   PetscScalar    v[3], stwo = -2./(appctx->h*appctx->h), sone = -.5*stwo;

469:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
470:      Compute entries for the locally owned part of the matrix
471:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
472:   /* 
473:      Set matrix rows corresponding to boundary data
474:   */

476:   mstart = 0;
477:   v[0] = 1.0;
478:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
479:   mstart++;

481:   mend--;
482:   v[0] = 1.0;
483:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

485:   /*
486:      Set matrix rows corresponding to interior data.  We construct the 
487:      matrix one row at a time.
488:   */
489:   v[0] = sone; v[1] = stwo; v[2] = sone;
490:   for ( i=mstart; i<mend; i++ ) {
491:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
492:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
493:   }

495:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
496:      Complete the matrix assembly process and set some options
497:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
498:   /*
499:      Assemble matrix, using the 2-step process:
500:        MatAssemblyBegin(), MatAssemblyEnd()
501:      Computations can be done while messages are in transition
502:      by placing code between these two statements.
503:   */
504:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
505:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

526:   /*
527:      Set and option to indicate that we will never add a new nonzero location 
528:      to the matrix. If we do, it will generate an error.
529:   */
530:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);

532:   return 0;
533: }
534: /* --------------------------------------------------------------------- */
537: /*
538:    Input Parameters:
539:    ts - the TS context
540:    t - current time
541:    f - function
542:    ctx - optional user-defined context, as set by TSetBCFunction()
543:  */
544: PetscErrorCode MyBCRoutine(TS ts,PetscReal t,Vec f,void *ctx)
545: {
546:   AppCtx         *appctx = (AppCtx *) ctx;     /* user-defined application context */
548:   PetscInt       m = appctx->m;
549:   PetscScalar    *fa;

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