Actual source code: ex5.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:   -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) = 1, u(t,1) = 1,
 25:    and the initial condition
 26:        u(0,x) = cos(6*pi*x) + 3*cos(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) * cos(6*pi*x) +
 38:                       3*exp(-4*pi*pi*t) * cos(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 just 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:      petsc.h       - base PETSc routines   petscvec.h  - vectors
 54:      petscsys.h    - system routines       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:   PetscTruth  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: */

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

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

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

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

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

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

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

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

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

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

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

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

156:   MatCreate(PETSC_COMM_SELF,&A);
157:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
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:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181:      Set solution vector and initial timestep
182:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

198:   TSSetDuration(ts,time_steps_max,time_total_max);
199:   TSSetFromOptions(ts);

201:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202:      Solve the problem
203:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

215:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216:      View timestepping solver info
217:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

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

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

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

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

283:   /* 
284:      Restore vector
285:   */
286:   VecRestoreArray(u,&u_localptr);

288:   /* 
289:      Print debugging information if desired
290:   */
291:   if (appctx->debug) {
292:      printf("initial guess vector\n");
293:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
294:   }

296:   return 0;
297: }
298: /* --------------------------------------------------------------------- */
301: /*
302:    ExactSolution - Computes the exact solution at a given time.

304:    Input Parameters:
305:    t - current time
306:    solution - vector in which exact solution will be computed
307:    appctx - user-defined application context

309:    Output Parameter:
310:    solution - vector with the newly computed exact solution
311: */
312: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
313: {
314:   PetscScalar    *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
315:   PetscInt       i;

318:   /*
319:      Get a pointer to vector data.
320:   */
321:   VecGetArray(solution,&s_localptr);

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

333:   /* 
334:      Restore vector
335:   */
336:   VecRestoreArray(solution,&s_localptr);
337:   return 0;
338: }
339: /* --------------------------------------------------------------------- */
342: /*
343:    Monitor - User-provided routine to monitor the solution computed at 
344:    each timestep.  This example plots the solution and computes the
345:    error in two different norms.

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

364:   /* 
365:      View a graph of the current iterate
366:   */
367:   VecView(u,appctx->viewer2);

369:   /* 
370:      Compute the exact solution
371:   */
372:   ExactSolution(time,appctx->solution,appctx);

374:   /*
375:      Print debugging information if desired
376:   */
377:   if (appctx->debug) {
378:      printf("Computed solution vector\n");
379:      VecView(u,PETSC_VIEWER_STDOUT_SELF);
380:      printf("Exact solution vector\n");
381:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
382:   }

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

392:   PetscPrintf(PETSC_COMM_WORLD,"Timestep %D: time = %G, 2-norm error = %G, max norm error = %G\n",
393:          step,time,norm_2,norm_max);
394:   appctx->norm_2   += norm_2;
395:   appctx->norm_max += norm_max;

397:   /* 
398:      View a graph of the error
399:   */
400:   VecView(appctx->solution,appctx->viewer1);

402:   /*
403:      Print debugging information if desired
404:   */
405:   if (appctx->debug) {
406:      printf("Error vector\n");
407:      VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
408:   }

410:   return 0;
411: }
412: /* --------------------------------------------------------------------- */
415: /*
416:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
417:    matrix for the heat equation.

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

425:    Output Parameters:
426:    AA - Jacobian matrix
427:    BB - optionally different preconditioning matrix
428:    str - flag indicating matrix structure

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

444:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
445:      Compute entries for the locally owned part of the matrix
446:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
447:   /* 
448:      Set matrix rows corresponding to boundary data
449:   */

451:   mstart = 0;
452:   v[0] = 1.0;
453:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
454:   mstart++;

456:   mend--;
457:   v[0] = 1.0;
458:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

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

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

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

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

507:   return 0;
508: }