Actual source code: ex2.c

petsc-master 2017-02-23
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*);

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

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:      Initialize program and set problem parameters
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 90:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 91:   PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,bounds);

 93:   appctx.comm = PETSC_COMM_WORLD;
 94:   appctx.m    = 60;

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

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

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:      Create vector data structures
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

106:   /*
107:      Create distributed array (DMDA) to manage parallel grid and vectors
108:      and to set up the ghost point communication pattern.  There are M
109:      total grid values spread equally among all the processors.
110:   */
111:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da);
112:   DMSetFromOptions(appctx.da);
113:   DMSetUp(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 ierr;
217: }
218: /* --------------------------------------------------------------------- */
219: /*
220:    InitialConditions - Computes the solution at the initial time.

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

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

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

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

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

262:   /*
263:      Restore vector
264:   */
265:   VecRestoreArray(u,&u_localptr);

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

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

281:    Input Parameters:
282:    t - current time
283:    solution - vector in which exact solution will be computed
284:    appctx - user-defined application context

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

295:   /*
296:      Determine starting and ending points of each processor's
297:      range of grid values
298:   */
299:   VecGetOwnershipRange(solution,&mybase,&myend);

301:   /*
302:      Get a pointer to vector data.
303:   */
304:   VecGetArray(solution,&s_localptr);

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

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

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

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

354:      PetscReal buffering makes graphics look better.
355:   */
356:   PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm),0,&draw);
357:   PetscDrawSetDoubleBuffer(draw);
358:   VecView(u,PETSC_VIEWER_DRAW_(appctx->comm));

360:   /*
361:      Compute the exact solution at this timestep
362:   */
363:   ExactSolution(time,appctx->solution,appctx);

365:   /*
366:      Print debugging information if desired
367:   */
368:   if (appctx->debug) {
369:     PetscPrintf(appctx->comm,"Computed solution vector\n");
370:     VecView(u,PETSC_VIEWER_STDOUT_WORLD);
371:     PetscPrintf(appctx->comm,"Exact solution vector\n");
372:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);
373:   }

375:   /*
376:      Compute the 2-norm and max-norm of the error
377:   */
378:   VecAXPY(appctx->solution,-1.0,u);
379:   VecNorm(appctx->solution,NORM_2,&en2);
380:   en2s = PetscSqrtReal(appctx->h)*en2;  /* scale the 2-norm by the grid spacing */
381:   VecNorm(appctx->solution,NORM_MAX,&enmax);

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

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

405:    Input Parameters:
406:    ts         - timesteping context
407:    t          - current time
408:    global_in  - vector containing the current iterate
409:    ctx        - (optional) user-provided context for function evaluation.
410:                 In this case we use the appctx defined above.

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

427:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
428:      Get ready for local function computations
429:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
430:   /*
431:      Scatter ghost points to local vector, using the 2-step process
432:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
433:      By placing code between these two statements, computations can be
434:      done while messages are in transition.
435:   */
436:   DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);
437:   DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);

439:   /*
440:       Access directly the values in our local INPUT work array
441:   */
442:   VecGetArrayRead(local_in,&localptr);

444:   /*
445:       Access directly the values in our local OUTPUT work array
446:   */
447:   VecGetArray(localwork,&copyptr);

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

451:   /*
452:       Evaluate our function on the nodes owned by this processor
453:   */
454:   VecGetLocalSize(local_in,&localsize);

456:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
457:      Compute entries for the locally owned part
458:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

460:   /*
461:      Handle boundary conditions: This is done by using the boundary condition
462:         u(t,boundary) = g(t,boundary)
463:      for some function g. Now take the derivative with respect to t to obtain
464:         u_{t}(t,boundary) = g_{t}(t,boundary)

466:      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
467:              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
468:   */
469:   MPI_Comm_rank(appctx->comm,&rank);
470:   MPI_Comm_size(appctx->comm,&size);
471:   if (!rank)          copyptr[0]           = 1.0;
472:   if (rank == size-1) copyptr[localsize-1] = 2.0;

474:   /*
475:      Handle the interior nodes where the PDE is replace by finite
476:      difference operators.
477:   */
478:   for (i=1; i<localsize-1; i++) copyptr[i] =  localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);

480:   /*
481:      Restore vectors
482:   */
483:   VecRestoreArrayRead(local_in,&localptr);
484:   VecRestoreArray(localwork,&copyptr);

486:   /*
487:      Insert values from the local OUTPUT vector into the global
488:      output vector
489:   */
490:   DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out);
491:   DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out);

493:   /* Print debugging information if desired */
494:   if (appctx->debug) {
495:     PetscPrintf(appctx->comm,"RHS function vector\n");
496:     VecView(global_out,PETSC_VIEWER_STDOUT_WORLD);
497:   }

499:   return 0;
500: }
501: /* --------------------------------------------------------------------- */
502: /*
503:    RHSJacobian - User-provided routine to compute the Jacobian of
504:    the nonlinear right-hand-side function of the ODE.

506:    Input Parameters:
507:    ts - the TS context
508:    t - current time
509:    global_in - global input vector
510:    dummy - optional user-defined context, as set by TSetRHSJacobian()

512:    Output Parameters:
513:    AA - Jacobian matrix
514:    BB - optionally different preconditioning matrix
515:    str - flag indicating matrix structure

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

540:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
541:      Get ready for local Jacobian computations
542:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
543:   /*
544:      Scatter ghost points to local vector, using the 2-step process
545:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
546:      By placing code between these two statements, computations can be
547:      done while messages are in transition.
548:   */
549:   DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);
550:   DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);

552:   /*
553:      Get pointer to vector data
554:   */
555:   VecGetArrayRead(local_in,&localptr);

557:   /*
558:      Get starting and ending locally owned rows of the matrix
559:   */
560:   MatGetOwnershipRange(BB,&mstarts,&mends);
561:   mstart = mstarts; mend = mends;

563:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
564:      Compute entries for the locally owned part of the Jacobian.
565:       - Currently, all PETSc parallel matrix formats are partitioned by
566:         contiguous chunks of rows across the processors.
567:       - Each processor needs to insert only elements that it owns
568:         locally (but any non-local elements will be sent to the
569:         appropriate processor during matrix assembly).
570:       - Here, we set all entries for a particular row at once.
571:       - We can set matrix entries either using either
572:         MatSetValuesLocal() or MatSetValues().
573:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

575:   /*
576:      Set matrix rows corresponding to boundary data
577:   */
578:   if (mstart == 0) {
579:     v[0] = 0.0;
580:     MatSetValues(BB,1,&mstart,1,&mstart,v,INSERT_VALUES);
581:     mstart++;
582:   }
583:   if (mend == appctx->m) {
584:     mend--;
585:     v[0] = 0.0;
586:     MatSetValues(BB,1,&mend,1,&mend,v,INSERT_VALUES);
587:   }

589:   /*
590:      Set matrix rows corresponding to interior data.  We construct the
591:      matrix one row at a time.
592:   */
593:   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
594:   for (i=mstart; i<mend; i++) {
595:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
596:     is     = i - mstart + 1;
597:     v[0]   = sc*localptr[is];
598:     v[1]   = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]);
599:     v[2]   = sc*localptr[is];
600:     MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES);
601:   }

603:   /*
604:      Restore vector
605:   */
606:   VecRestoreArrayRead(local_in,&localptr);

608:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
609:      Complete the matrix assembly process and set some options
610:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
611:   /*
612:      Assemble matrix, using the 2-step process:
613:        MatAssemblyBegin(), MatAssemblyEnd()
614:      Computations can be done while messages are in transition
615:      by placing code between these two statements.
616:   */
617:   MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);
618:   MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);
619:   if (BB != AA) {
620:     MatAssemblyBegin(AA,MAT_FINAL_ASSEMBLY);
621:     MatAssemblyEnd(AA,MAT_FINAL_ASSEMBLY);
622:   }

624:   /*
625:      Set and option to indicate that we will never add a new nonzero location
626:      to the matrix. If we do, it will generate an error.
627:   */
628:   MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

630:   return 0;
631: }