Actual source code: ex2.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] ="Solves a time-dependent nonlinear PDE. Uses implicit\n\
  3: timestepping.  Runtime options include:\n\
  4:   -M <xg>, where <xg> = number of grid points\n\
  5:   -debug : Activate debugging printouts\n\
  6:   -nox   : Deactivate x-window graphics\n\n";

  8: /*
  9:    Concepts: TS^time-dependent nonlinear problems
 10:    Processors: n
 11: */

 13: /* ------------------------------------------------------------------------

 15:    This program solves the PDE

 17:                u * u_xx
 18:          u_t = ---------
 19:                2*(t+1)^2

 21:     on the domain 0 <= x <= 1, with boundary conditions
 22:          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
 23:     and initial condition
 24:          u(0,x) = 1 + x*x.

 26:     The exact solution is:
 27:          u(t,x) = (1 + x*x) * (1 + t)

 29:     Note that since the solution is linear in time and quadratic in x,
 30:     the finite difference scheme actually computes the "exact" solution.

 32:     We use by default the backward Euler method.

 34:   ------------------------------------------------------------------------- */

 36: /*
 37:    Include "petscts.h" to use the PETSc timestepping routines. Note that
 38:    this file automatically includes "petscsys.h" and other lower-level
 39:    PETSc include files.

 41:    Include the "petscdmda.h" to allow us to use the distributed array data
 42:    structures to manage the parallel grid.
 43: */
 44: #include <petscts.h>
 45: #include <petscdm.h>
 46: #include <petscdmda.h>
 47: #include <petscdraw.h>

 49: /*
 50:    User-defined application context - contains data needed by the
 51:    application-provided callback routines.
 52: */
 53: typedef struct {
 54:   MPI_Comm  comm;           /* communicator */
 55:   DM        da;             /* distributed array data structure */
 56:   Vec       localwork;      /* local ghosted work vector */
 57:   Vec       u_local;        /* local ghosted approximate solution vector */
 58:   Vec       solution;       /* global exact solution vector */
 59:   PetscInt  m;              /* total number of grid points */
 60:   PetscReal h;              /* mesh width: h = 1/(m-1) */
 61:   PetscBool debug;          /* flag (1 indicates activation of debugging printouts) */
 62: } AppCtx;

 64: /*
 65:    User-defined routines, provided below.
 66: */
 67: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
 68: extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
 69: extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
 70: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 71: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);

 75: int main(int argc,char **argv)
 76: {
 77:   AppCtx         appctx;                 /* user-defined application context */
 78:   TS             ts;                     /* timestepping context */
 79:   Mat            A;                      /* Jacobian matrix data structure */
 80:   Vec            u;                      /* approximate solution vector */
 81:   PetscInt       time_steps_max = 100;  /* default max timesteps */
 83:   PetscReal      dt;
 84:   PetscReal      time_total_max = 100.0; /* default max total time */
 85:   PetscBool      mymonitor      = PETSC_FALSE;
 86:   PetscReal      bounds[]       = {1.0, 3.3};

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:      Initialize program and set problem parameters
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   PetscInitialize(&argc,&argv,(char*)0,help);
 93:   PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,bounds);

 95:   appctx.comm = PETSC_COMM_WORLD;
 96:   appctx.m    = 60;

 98:   PetscOptionsGetInt(NULL,NULL,"-M",&appctx.m,NULL);
 99:   PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);
100:   PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);

102:   appctx.h    = 1.0/(appctx.m-1.0);

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:      Create vector data structures
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   /*
109:      Create distributed array (DMDA) to manage parallel grid and vectors
110:      and to set up the ghost point communication pattern.  There are M
111:      total grid values spread equally among all the processors.
112:   */
113:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da);

115:   /*
116:      Extract global and local vectors from DMDA; we use these to store the
117:      approximate solution.  Then duplicate these for remaining vectors that
118:      have the same types.
119:   */
120:   DMCreateGlobalVector(appctx.da,&u);
121:   DMCreateLocalVector(appctx.da,&appctx.u_local);

123:   /*
124:      Create local work vector for use in evaluating right-hand-side function;
125:      create global work vector for storing exact solution.
126:   */
127:   VecDuplicate(appctx.u_local,&appctx.localwork);
128:   VecDuplicate(u,&appctx.solution);

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:      Create timestepping solver context; set callback routine for
132:      right-hand-side function evaluation.
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

135:   TSCreate(PETSC_COMM_WORLD,&ts);
136:   TSSetProblemType(ts,TS_NONLINEAR);
137:   TSSetRHSFunction(ts,NULL,RHSFunction,&appctx);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:      Set optional user-defined monitoring routine
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

143:   if (mymonitor) {
144:     TSMonitorSet(ts,Monitor,&appctx,NULL);
145:   }

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:      For nonlinear problems, the user can provide a Jacobian evaluation
149:      routine (or use a finite differencing approximation).

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

154:   MatCreate(PETSC_COMM_WORLD,&A);
155:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m);
156:   MatSetFromOptions(A);
157:   MatSetUp(A);
158:   TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx);

160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161:      Set solution vector and initial timestep
162:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

164:   dt   = appctx.h/2.0;
165:   TSSetInitialTimeStep(ts,0.0,dt);

167:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168:      Customize timestepping solver:
169:        - Set the solution method to be the Backward Euler method.
170:        - Set timestepping duration info
171:      Then set runtime options, which can override these defaults.
172:      For example,
173:           -ts_max_steps <maxsteps> -ts_final_time <maxtime>
174:      to override the defaults set by TSSetDuration().
175:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

177:   TSSetType(ts,TSBEULER);
178:   TSSetDuration(ts,time_steps_max,time_total_max);
179:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
180:   TSSetFromOptions(ts);
181: 
182:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183:      Solve the problem
184:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

186:   /*
187:      Evaluate initial conditions
188:   */
189:   InitialConditions(u,&appctx);

191:   /*
192:      Run the timestepping solver
193:   */
194:   TSSolve(ts,u);

196:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197:      Free work space.  All PETSc objects should be destroyed when they
198:      are no longer needed.
199:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

201:   TSDestroy(&ts);
202:   VecDestroy(&u);
203:   MatDestroy(&A);
204:   DMDestroy(&appctx.da);
205:   VecDestroy(&appctx.localwork);
206:   VecDestroy(&appctx.solution);
207:   VecDestroy(&appctx.u_local);

209:   /*
210:      Always call PetscFinalize() before exiting a program.  This routine
211:        - finalizes the PETSc libraries as well as MPI
212:        - provides summary and diagnostic information if certain runtime
213:          options are chosen (e.g., -log_summary).
214:   */
215:   PetscFinalize();
216:   return 0;
217: }
218: /* --------------------------------------------------------------------- */
221: /*
222:    InitialConditions - Computes the solution at the initial time.

224:    Input Parameters:
225:    u - uninitialized solution vector (global)
226:    appctx - user-defined application context

228:    Output Parameter:
229:    u - vector with solution at initial time (global)
230: */
231: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
232: {
233:   PetscScalar    *u_localptr,h = appctx->h,x;
234:   PetscInt       i,mybase,myend;

237:   /*
238:      Determine starting point of each processor's range of
239:      grid values.
240:   */
241:   VecGetOwnershipRange(u,&mybase,&myend);

243:   /*
244:     Get a pointer to vector data.
245:     - For default PETSc vectors, VecGetArray() returns a pointer to
246:       the data array.  Otherwise, the routine is implementation dependent.
247:     - You MUST call VecRestoreArray() when you no longer need access to
248:       the array.
249:     - Note that the Fortran interface to VecGetArray() differs from the
250:       C version.  See the users manual for details.
251:   */
252:   VecGetArray(u,&u_localptr);

254:   /*
255:      We initialize the solution array by simply writing the solution
256:      directly into the array locations.  Alternatively, we could use
257:      VecSetValues() or VecSetValuesLocal().
258:   */
259:   for (i=mybase; i<myend; i++) {
260:     x = h*(PetscReal)i; /* current location in global grid */
261:     u_localptr[i-mybase] = 1.0 + x*x;
262:   }

264:   /*
265:      Restore vector
266:   */
267:   VecRestoreArray(u,&u_localptr);

269:   /*
270:      Print debugging information if desired
271:   */
272:   if (appctx->debug) {
273:     PetscPrintf(appctx->comm,"initial guess vector\n");
274:     VecView(u,PETSC_VIEWER_STDOUT_WORLD);
275:   }

277:   return 0;
278: }
279: /* --------------------------------------------------------------------- */
282: /*
283:    ExactSolution - Computes the exact solution at a given time.

285:    Input Parameters:
286:    t - current time
287:    solution - vector in which exact solution will be computed
288:    appctx - user-defined application context

290:    Output Parameter:
291:    solution - vector with the newly computed exact solution
292: */
293: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
294: {
295:   PetscScalar    *s_localptr,h = appctx->h,x;
296:   PetscInt       i,mybase,myend;

299:   /*
300:      Determine starting and ending points of each processor's
301:      range of grid values
302:   */
303:   VecGetOwnershipRange(solution,&mybase,&myend);

305:   /*
306:      Get a pointer to vector data.
307:   */
308:   VecGetArray(solution,&s_localptr);

310:   /*
311:      Simply write the solution directly into the array locations.
312:      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
313:   */
314:   for (i=mybase; i<myend; i++) {
315:     x = h*(PetscReal)i;
316:     s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x);
317:   }

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

333:    Input Parameters:
334:    ts     - the timestep context
335:    step   - the count of the current step (with 0 meaning the
336:             initial condition)
337:    time   - the current time
338:    u      - the solution at this timestep
339:    ctx    - the user-provided context for this monitoring routine.
340:             In this case we use the application context which contains
341:             information about the problem size, workspace and the exact
342:             solution.
343: */
344: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
345: {
346:   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
348:   PetscReal      en2,en2s,enmax;
349:   PetscDraw      draw;

351:   /*
352:      We use the default X windows viewer
353:              PETSC_VIEWER_DRAW_(appctx->comm)
354:      that is associated with the current communicator. This saves
355:      the effort of calling PetscViewerDrawOpen() to create the window.
356:      Note that if we wished to plot several items in separate windows we
357:      would create each viewer with PetscViewerDrawOpen() and store them in
358:      the application context, appctx.

360:      PetscReal buffering makes graphics look better.
361:   */
362:   PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm),0,&draw);
363:   PetscDrawSetDoubleBuffer(draw);
364:   VecView(u,PETSC_VIEWER_DRAW_(appctx->comm));

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

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

381:   /*
382:      Compute the 2-norm and max-norm of the error
383:   */
384:   VecAXPY(appctx->solution,-1.0,u);
385:   VecNorm(appctx->solution,NORM_2,&en2);
386:   en2s = PetscSqrtReal(appctx->h)*en2;  /* scale the 2-norm by the grid spacing */
387:   VecNorm(appctx->solution,NORM_MAX,&enmax);

389:   /*
390:      PetscPrintf() causes only the first processor in this
391:      communicator to print the timestep information.
392:   */
393:   PetscPrintf(appctx->comm,"Timestep %D: time = %g 2-norm error = %g  max norm error = %g\n",step,(double)time,(double)en2s,(double)enmax);

395:   /*
396:      Print debugging information if desired
397:   */
398:   if (appctx->debug) {
399:     PetscPrintf(appctx->comm,"Error vector\n");
400:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);
401:   }
402:   return 0;
403: }
404: /* --------------------------------------------------------------------- */
407: /*
408:    RHSFunction - User-provided routine that evalues the right-hand-side
409:    function of the ODE.  This routine is set in the main program by
410:    calling TSSetRHSFunction().  We compute:
411:           global_out = F(global_in)

413:    Input Parameters:
414:    ts         - timesteping context
415:    t          - current time
416:    global_in  - vector containing the current iterate
417:    ctx        - (optional) user-provided context for function evaluation.
418:                 In this case we use the appctx defined above.

420:    Output Parameter:
421:    global_out - vector containing the newly evaluated function
422: */
423: PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx)
424: {
425:   AppCtx            *appctx   = (AppCtx*) ctx;     /* user-defined application context */
426:   DM                da        = appctx->da;        /* distributed array */
427:   Vec               local_in  = appctx->u_local;   /* local ghosted input vector */
428:   Vec               localwork = appctx->localwork; /* local ghosted work vector */
429:   PetscErrorCode    ierr;
430:   PetscInt          i,localsize;
431:   PetscMPIInt       rank,size;
432:   PetscScalar       *copyptr,sc;
433:   const PetscScalar *localptr;

435:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
436:      Get ready for local function computations
437:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
438:   /*
439:      Scatter ghost points to local vector, using the 2-step process
440:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
441:      By placing code between these two statements, computations can be
442:      done while messages are in transition.
443:   */
444:   DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);
445:   DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);

447:   /*
448:       Access directly the values in our local INPUT work array
449:   */
450:   VecGetArrayRead(local_in,&localptr);

452:   /*
453:       Access directly the values in our local OUTPUT work array
454:   */
455:   VecGetArray(localwork,&copyptr);

457:   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));

459:   /*
460:       Evaluate our function on the nodes owned by this processor
461:   */
462:   VecGetLocalSize(local_in,&localsize);

464:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
465:      Compute entries for the locally owned part
466:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

468:   /*
469:      Handle boundary conditions: This is done by using the boundary condition
470:         u(t,boundary) = g(t,boundary)
471:      for some function g. Now take the derivative with respect to t to obtain
472:         u_{t}(t,boundary) = g_{t}(t,boundary)

474:      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
475:              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
476:   */
477:   MPI_Comm_rank(appctx->comm,&rank);
478:   MPI_Comm_size(appctx->comm,&size);
479:   if (!rank)          copyptr[0]           = 1.0;
480:   if (rank == size-1) copyptr[localsize-1] = 2.0;

482:   /*
483:      Handle the interior nodes where the PDE is replace by finite
484:      difference operators.
485:   */
486:   for (i=1; i<localsize-1; i++) copyptr[i] =  localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);

488:   /*
489:      Restore vectors
490:   */
491:   VecRestoreArrayRead(local_in,&localptr);
492:   VecRestoreArray(localwork,&copyptr);

494:   /*
495:      Insert values from the local OUTPUT vector into the global
496:      output vector
497:   */
498:   DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out);
499:   DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out);

501:   /* Print debugging information if desired */
502:   if (appctx->debug) {
503:     PetscPrintf(appctx->comm,"RHS function vector\n");
504:     VecView(global_out,PETSC_VIEWER_STDOUT_WORLD);
505:   }

507:   return 0;
508: }
509: /* --------------------------------------------------------------------- */
512: /*
513:    RHSJacobian - User-provided routine to compute the Jacobian of
514:    the nonlinear right-hand-side function of the ODE.

516:    Input Parameters:
517:    ts - the TS context
518:    t - current time
519:    global_in - global input vector
520:    dummy - optional user-defined context, as set by TSetRHSJacobian()

522:    Output Parameters:
523:    AA - Jacobian matrix
524:    BB - optionally different preconditioning matrix
525:    str - flag indicating matrix structure

527:   Notes:
528:   RHSJacobian computes entries for the locally owned part of the Jacobian.
529:    - Currently, all PETSc parallel matrix formats are partitioned by
530:      contiguous chunks of rows across the processors.
531:    - Each processor needs to insert only elements that it owns
532:      locally (but any non-local elements will be sent to the
533:      appropriate processor during matrix assembly).
534:    - Always specify global row and columns of matrix entries when
535:      using MatSetValues().
536:    - Here, we set all entries for a particular row at once.
537:    - Note that MatSetValues() uses 0-based row and column numbers
538:      in Fortran as well as in C.
539: */
540: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat BB,void *ctx)
541: {
542:   AppCtx            *appctx  = (AppCtx*)ctx;    /* user-defined application context */
543:   Vec               local_in = appctx->u_local;   /* local ghosted input vector */
544:   DM                da       = appctx->da;        /* distributed array */
545:   PetscScalar       v[3],sc;
546:   const PetscScalar *localptr;
547:   PetscErrorCode    ierr;
548:   PetscInt          i,mstart,mend,mstarts,mends,idx[3],is;

550:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
551:      Get ready for local Jacobian computations
552:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
553:   /*
554:      Scatter ghost points to local vector, using the 2-step process
555:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
556:      By placing code between these two statements, computations can be
557:      done while messages are in transition.
558:   */
559:   DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);
560:   DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);

562:   /*
563:      Get pointer to vector data
564:   */
565:   VecGetArrayRead(local_in,&localptr);

567:   /*
568:      Get starting and ending locally owned rows of the matrix
569:   */
570:   MatGetOwnershipRange(BB,&mstarts,&mends);
571:   mstart = mstarts; mend = mends;

573:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
574:      Compute entries for the locally owned part of the Jacobian.
575:       - Currently, all PETSc parallel matrix formats are partitioned by
576:         contiguous chunks of rows across the processors.
577:       - Each processor needs to insert only elements that it owns
578:         locally (but any non-local elements will be sent to the
579:         appropriate processor during matrix assembly).
580:       - Here, we set all entries for a particular row at once.
581:       - We can set matrix entries either using either
582:         MatSetValuesLocal() or MatSetValues().
583:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

585:   /*
586:      Set matrix rows corresponding to boundary data
587:   */
588:   if (mstart == 0) {
589:     v[0] = 0.0;
590:     MatSetValues(BB,1,&mstart,1,&mstart,v,INSERT_VALUES);
591:     mstart++;
592:   }
593:   if (mend == appctx->m) {
594:     mend--;
595:     v[0] = 0.0;
596:     MatSetValues(BB,1,&mend,1,&mend,v,INSERT_VALUES);
597:   }

599:   /*
600:      Set matrix rows corresponding to interior data.  We construct the
601:      matrix one row at a time.
602:   */
603:   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
604:   for (i=mstart; i<mend; i++) {
605:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
606:     is     = i - mstart + 1;
607:     v[0]   = sc*localptr[is];
608:     v[1]   = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]);
609:     v[2]   = sc*localptr[is];
610:     MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES);
611:   }

613:   /*
614:      Restore vector
615:   */
616:   VecRestoreArrayRead(local_in,&localptr);

618:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
619:      Complete the matrix assembly process and set some options
620:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
621:   /*
622:      Assemble matrix, using the 2-step process:
623:        MatAssemblyBegin(), MatAssemblyEnd()
624:      Computations can be done while messages are in transition
625:      by placing code between these two statements.
626:   */
627:   MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);
628:   MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);
629:   if (BB != AA) {
630:     MatAssemblyBegin(AA,MAT_FINAL_ASSEMBLY);
631:     MatAssemblyEnd(AA,MAT_FINAL_ASSEMBLY);
632:   }

634:   /*
635:      Set and option to indicate that we will never add a new nonzero location
636:      to the matrix. If we do, it will generate an error.
637:   */
638:   MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

640:   return 0;
641: }