Actual source code: ex3.c

petsc-master 2017-06-23
Report Typos and Errors

  2: static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\
  3: Input parameters include:\n\
  4:   -m <points>, where <points> = number of grid points\n\
  5:   -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
  6:   -debug              : Activate debugging printouts\n\
  7:   -nox                : Deactivate x-window graphics\n\n";

  9: /*
 10:    Concepts: TS^time-dependent linear problems
 11:    Concepts: TS^heat equation
 12:    Concepts: TS^diffusion equation
 13:    Processors: 1
 14: */

 16: /* ------------------------------------------------------------------------

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

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

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

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

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

 46:   ------------------------------------------------------------------------- */

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

 58:  #include <petscts.h>
 59:  #include <petscdraw.h>

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

 74: /*
 75:    User-defined routines
 76: */
 77: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
 78: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
 79: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 80: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);

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

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:      Initialize program and set problem parameters
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

101:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
102:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
103:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

105:   m    = 60;
106:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
107:   PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);

109:   appctx.m        = m;
110:   appctx.h        = 1.0/(m-1.0);
111:   appctx.norm_2   = 0.0;
112:   appctx.norm_max = 0.0;

114:   PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");

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

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

155:   MatCreate(PETSC_COMM_SELF,&A);
156:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
157:   MatSetFromOptions(A);
158:   MatSetUp(A);

160:   flg  = PETSC_FALSE;
161:   PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL);
162:   if (flg) {
163:     /*
164:        For linear problems with a time-dependent f(u,t) in the equation
165:        u_t = f(u,t), the user provides the discretized right-hand-side
166:        as a time-dependent matrix.
167:     */
168:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
169:     TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
170:   } else {
171:     /*
172:        For linear problems with a time-independent f(u) in the equation
173:        u_t = f(u), the user provides the discretized right-hand-side
174:        as a matrix only once, and then sets the special Jacobian evaluation
175:        routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian.
176:     */
177:     RHSMatrixHeat(ts,0.0,u,A,A,&appctx);
178:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
179:     TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
180:   }

182:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183:      Set solution vector and initial timestep
184:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

186:   dt   = appctx.h*appctx.h/2.0;
187:   TSSetInitialTimeStep(ts,0.0,dt);

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

199:   TSSetDuration(ts,time_steps_max,time_total_max);
200:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
201:   TSSetFromOptions(ts);

203:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204:      Solve the problem
205:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

212:   /*
213:      Run the timestepping solver
214:   */
215:   TSSolve(ts,u);
216:   TSGetTimeStepNumber(ts,&steps);

218:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219:      View timestepping solver info
220:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

222:   PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %g, avg. error (max norm) = %g\n",(double)(appctx.norm_2/steps),(double)(appctx.norm_max/steps));
223:   TSView(ts,PETSC_VIEWER_STDOUT_SELF);

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

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

237:   /*
238:      Always call PetscFinalize() before exiting a program.  This routine
239:        - finalizes the PETSc libraries as well as MPI
240:        - provides summary and diagnostic information if certain runtime
241:          options are chosen (e.g., -log_view).
242:   */
243:   PetscFinalize();
244:   return ierr;
245: }
246: /* --------------------------------------------------------------------- */
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;
261:   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++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);

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

286:   /*
287:      Print debugging information if desired
288:   */
289:   if (appctx->debug) {
290:     PetscPrintf(PETSC_COMM_WORLD,"Initial guess vector\n");
291:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
292:   }

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

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

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

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);
324:   ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
325:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
326:   for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;

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

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

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

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

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

369:   /*
370:      Print debugging information if desired
371:   */
372:   if (appctx->debug) {
373:     PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n");
374:     VecView(u,PETSC_VIEWER_STDOUT_SELF);
375:     PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");
376:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
377:   }

379:   /*
380:      Compute the 2-norm and max-norm of the error
381:   */
382:   VecAXPY(appctx->solution,-1.0,u);
383:   VecNorm(appctx->solution,NORM_2,&norm_2);
384:   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
385:   VecNorm(appctx->solution,NORM_MAX,&norm_max);
386:   if (norm_2   < 1e-14) norm_2   = 0;
387:   if (norm_max < 1e-14) norm_max = 0;

389:   TSGetTimeStep(ts,&dt);
390:   PetscPrintf(PETSC_COMM_WORLD,"Timestep %3D: step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)dt,(double)time,(double)norm_2,(double)norm_max);

392:   appctx->norm_2   += norm_2;
393:   appctx->norm_max += norm_max;

395:   dttol = .0001;
396:   PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,NULL);
397:   if (dt < dttol) {
398:     dt  *= .999;
399:     TSSetTimeStep(ts,dt);
400:   }

402:   /*
403:      View a graph of the error
404:   */
405:   VecView(appctx->solution,appctx->viewer1);

407:   /*
408:      Print debugging information if desired
409:   */
410:   if (appctx->debug) {
411:     PetscPrintf(PETSC_COMM_SELF,"Error vector\n");
412:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);
413:   }

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

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

428:    Output Parameters:
429:    AA - Jacobian matrix
430:    BB - optionally different preconditioning matrix
431:    str - flag indicating matrix structure

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

447:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
448:      Compute entries for the locally owned part of the matrix
449:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
450:   /*
451:      Set matrix rows corresponding to boundary data
452:   */

454:   mstart = 0;
455:   v[0]   = 1.0;
456:   MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
457:   mstart++;

459:   mend--;
460:   v[0] = 1.0;
461:   MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);

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

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

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

491:   return 0;
492: }