Actual source code: ex14f.F

  1: !
  2: !
  3: !  Solves a nonlinear system in parallel with a user-defined
  4: !  Newton method that uses KSP to solve the linearized Newton sytems.  This solver
  5: !  is a very simplistic inexact Newton method.  The intent of this code is to
  6: !  demonstrate the repeated solution of linear sytems with the same nonzero pattern.
  7: !
  8: !  This is NOT the recommended approach for solving nonlinear problems with PETSc!
  9: !  We urge users to employ the SNES component for solving nonlinear problems whenever
 10: !  possible, as it offers many advantages over coding nonlinear solvers independently.
 11: !
 12: !  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
 13: !  domain, using distributed arrays (DAs) to partition the parallel grid.
 14: !
 15: !  The command line options include:
 16: !  -par <parameter>, where <parameter> indicates the problem's nonlinearity
 17: !     problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
 18: !  -mx <xg>, where <xg> = number of grid points in the x-direction
 19: !  -my <yg>, where <yg> = number of grid points in the y-direction
 20: !  -Nx <npx>, where <npx> = number of processors in the x-direction
 21: !  -Ny <npy>, where <npy> = number of processors in the y-direction
 22: !  -mf use matrix free for matrix vector product
 23: !
 24: !/*T
 25: !   Concepts: KSP^writing a user-defined nonlinear solver
 26: !   Concepts: DA^using distributed arrays
 27: !   Processors: n
 28: !T*/
 29: !  ------------------------------------------------------------------------
 30: !
 31: !    Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 32: !    the partial differential equation
 33: !
 34: !            -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 35: !
 36: !    with boundary conditions
 37: !
 38: !             u = 0  for  x = 0, x = 1, y = 0, y = 1.
 39: !
 40: !    A finite difference approximation with the usual 5-point stencil
 41: !    is used to discretize the boundary value problem to obtain a nonlinear
 42: !    system of equations.
 43: !
 44: !    The SNES version of this problem is:  snes/examples/tutorials/ex5f.F
 45: !
 46: !  -------------------------------------------------------------------------

 48:       program main
 49:       implicit none

 51: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 52: !                    Include files
 53: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 54: !
 55: !     petsc.h       - base PETSc routines   petscvec.h - vectors
 56: !     petscsys.h    - system routines       petscmat.h - matrices
 57: !     petscis.h     - index sets            petscksp.h - Krylov subspace methods
 58: !     petscviewer.h - viewers               petscpc.h  - preconditioners

 60:  #include include/finclude/petsc.h
 61:  #include include/finclude/petscis.h
 62:  #include include/finclude/petscvec.h
 63:  #include include/finclude/petscmat.h
 64:  #include include/finclude/petscpc.h
 65:  #include include/finclude/petscksp.h
 66:  #include include/finclude/petscda.h

 68:       MPI_Comm comm
 69:       Vec      X,Y,F,localX,localF
 70:       Mat      J,B
 71:       DA       da
 72:       KSP      ksp

 74:       PetscInt  Nx,Ny,N,mx,my,ifive,ithree
 75:       PetscTruth flg,nooutput,usemf
 76:       common   /mycommon/ mx,my,B,localX,localF,da
 77: !
 78: !
 79: !      This is the routine to use for matrix-free approach
 80: !
 81:       external mymult

 83: !     --------------- Data to define nonlinear solver --------------
 84:       double precision   rtol,xtol,ttol
 85:       double precision   fnorm,ynorm,xnorm
 86:       PetscInt            max_nonlin_its,one
 87:       PetscInt            lin_its
 88:       PetscInt           i,m
 89:       PetscScalar        mone
 90:       PetscErrorCode ierr

 92:       mone           = -1.d0
 93:       rtol           = 1.d-8
 94:       xtol           = 1.d-8
 95:       max_nonlin_its = 10
 96:       one            = 1
 97:       ifive          = 5
 98:       ithree         = 3

100:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
101:       comm = PETSC_COMM_WORLD

103: !  Initialize problem parameters

105: !
106:       mx = 4
107:       my = 4
108:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr)
109:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-my',my,flg,ierr)
110:       N = mx*my

112:       nooutput = 0
113:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-no_output',       &
114:      &     nooutput,ierr)

116: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: !     Create linear solver context
118: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

120:       call KSPCreate(comm,ksp,ierr)

122: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: !     Create vector data structures
124: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

126: !
127: !  Create distributed array (DA) to manage parallel grid and vectors
128: !
129:       Nx = PETSC_DECIDE
130:       Ny = PETSC_DECIDE
131:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-Nx',Nx,flg,ierr)
132:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-Ny',Ny,flg,ierr)
133:       call DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_STAR,mx,           &
134:      &     my,Nx,Ny,one,one,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,       &
135:      &     da,ierr)

137: !
138: !  Extract global and local vectors from DA then duplicate for remaining
139: !  vectors that are the same types
140: !
141:        call DACreateGlobalVector(da,X,ierr)
142:        call DACreateLocalVector(da,localX,ierr)
143:        call VecDuplicate(X,F,ierr)
144:        call VecDuplicate(X,Y,ierr)
145:        call VecDuplicate(localX,localF,ierr)


148: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: !     Create matrix data structure for Jacobian
150: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: !
152: !     Note:  For the parallel case, vectors and matrices MUST be partitioned
153: !     accordingly.  When using distributed arrays (DAs) to create vectors,
154: !     the DAs determine the problem partitioning.  We must explicitly
155: !     specify the local matrix dimensions upon its creation for compatibility
156: !     with the vector distribution.
157: !
158: !     Note: Here we only approximately preallocate storage space for the
159: !     Jacobian.  See the users manual for a discussion of better techniques
160: !     for preallocating matrix memory.
161: !
162:       call VecGetLocalSize(X,m,ierr)
163:       call MatCreateMPIAIJ(comm,m,m,N,N,ifive,PETSC_NULL_INTEGER,ithree,         &
164:      &     PETSC_NULL_INTEGER,B,ierr)

166: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: !     if usemf is on then matrix vector product is done via matrix free
168: !     approach. Note this is just an example, and not realistic because
169: !     we still use the actual formed matrix, but in reality one would
170: !     provide their own subroutine that would directly do the matrix
171: !     vector product and not call MatMult()
172: !     Note: we put B into a common block so it will be visible to the
173: !     mymult() routine
174:       usemf = 0
175:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-mf',usemf,ierr)
176:       if (usemf .eq. 1) then
177:          call MatCreateShell(comm,m,m,N,N,PETSC_NULL_INTEGER,J,ierr)
178:          call MatShellSetOperation(J,MATOP_MULT,mymult,ierr)
179:       else
180: !        If not doing matrix free then matrix operator, J,  and matrix used
181: !        to construct preconditioner, B, are the same
182:         J = B
183:       endif

185: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186: !     Customize linear solver set runtime options
187: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188: !
189: !     Set runtime options (e.g., -ksp_monitor -ksp_rtol <rtol> -ksp_type <type>)
190: !
191:        call KSPSetFromOptions(ksp,ierr)

193: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194: !     Evaluate initial guess
195: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

197:        call FormInitialGuess(X,ierr)
198:        call ComputeFunction(X,F,ierr)
199:        call VecNorm(F,NORM_2,fnorm,ierr)
200:        ttol = fnorm*rtol
201:        if (nooutput .eq. 0) then
202:          print*, 'Initial function norm ',fnorm
203:        endif

205: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: !     Solve nonlinear system with a user-defined method
207: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

209: !  This solver is a very simplistic inexact Newton method, with no
210: !  no damping strategies or bells and whistles. The intent of this code
211: !  is merely to demonstrate the repeated solution with KSP of linear
212: !  sytems with the same nonzero structure.
213: !
214: !  This is NOT the recommended approach for solving nonlinear problems
215: !  with PETSc!  We urge users to employ the SNES component for solving
216: !  nonlinear problems whenever possible with application codes, as it
217: !  offers many advantages over coding nonlinear solvers independently.

219:        do 10 i=0,max_nonlin_its

221: !  Compute the Jacobian matrix.  See the comments in this routine for
222: !  important information about setting the flag mat_flag.

224:          call ComputeJacobian(X,B,ierr)

226: !  Solve J Y = F, where J is the Jacobian matrix.
227: !    - First, set the KSP linear operators.  Here the matrix that
228: !      defines the linear system also serves as the preconditioning
229: !      matrix.
230: !    - Then solve the Newton system.

232:          call KSPSetOperators(ksp,J,B,SAME_NONZERO_PATTERN,ierr)
233:          call KSPSolve(ksp,F,Y,ierr)

235: !  Compute updated iterate

237:          call VecNorm(Y,NORM_2,ynorm,ierr)
238:          call VecAYPX(Y,mone,X,ierr)
239:          call VecCopy(Y,X,ierr)
240:          call VecNorm(X,NORM_2,xnorm,ierr)
241:          call KSPGetIterationNumber(ksp,lin_its,ierr)
242:          if (nooutput .eq. 0) then
243:            print*,'linear solve iterations = ',lin_its,' xnorm = ',     &
244:      &         xnorm,' ynorm = ',ynorm
245:          endif

247: !  Evaluate nonlinear function at new location

249:          call ComputeFunction(X,F,ierr)
250:          call VecNorm(F,NORM_2,fnorm,ierr)
251:          if (nooutput .eq. 0) then
252:            print*, 'Iteration ',i+1,' function norm',fnorm
253:          endif

255: !  Test for convergence

257:        if (fnorm .le. ttol) then
258:          if (nooutput .eq. 0) then
259:            print*,'Converged: function norm ',fnorm,' tolerance ',ttol
260:          endif
261:          goto 20
262:        endif
263:  10   continue
264:  20   continue

266:       write(6,100) i+1
267:  100  format('Number of Newton iterations =',I2)

269: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270: !     Free work space.  All PETSc objects should be destroyed when they
271: !     are no longer needed.
272: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

274:        call MatDestroy(B,ierr)
275:        if (usemf .ne. 0) then
276:          call MatDestroy(J,ierr)
277:        endif
278:        call VecDestroy(localX,ierr)
279:        call VecDestroy(X,ierr)
280:        call VecDestroy(Y,ierr)
281:        call VecDestroy(localF,ierr)
282:        call VecDestroy(F,ierr)
283:        call KSPDestroy(ksp,ierr)
284:        call DADestroy(da,ierr)
285:        call PetscFinalize(ierr)
286:        end

288: ! -------------------------------------------------------------------
289: !
290: !   FormInitialGuess - Forms initial approximation.
291: !
292: !   Input Parameters:
293: !   X - vector
294: !
295: !   Output Parameter:
296: !   X - vector
297: !
298:       subroutine FormInitialGuess(X,ierr)
299:       implicit none

301: !     petsc.h       - base PETSc routines   petscvec.h - vectors
302: !     petscsys.h    - system routines       petscmat.h - matrices
303: !     petscis.h     - index sets            petscksp.h - Krylov subspace methods
304: !     petscviewer.h - viewers               petscpc.h  - preconditioners

306:  #include include/finclude/petsc.h
307:  #include include/finclude/petscis.h
308:  #include include/finclude/petscvec.h
309:  #include include/finclude/petscmat.h
310:  #include include/finclude/petscpc.h
311:  #include include/finclude/petscksp.h
312:  #include include/finclude/petscda.h
313:       PetscErrorCode    ierr
314:       PetscOffset      idx
315:       Vec       X,localX,localF
316:       PetscInt  i,j,row,mx
317:       PetscInt  my, xs,ys,xm
318:       PetscInt  ym,gxm,gym,gxs,gys
319:       double precision one,lambda,temp1,temp,hx,hy
320:       double precision hxdhy,hydhx,sc
321:       PetscScalar      xx(1)
322:       DA               da
323:       Mat              B
324:       common   /mycommon/ mx,my,B,localX,localF,da
325: 
326:       one    = 1.d0
327:       lambda = 6.d0
328:       hx     = one/(mx-1)
329:       hy     = one/(my-1)
330:       sc     = hx*hy*lambda
331:       hxdhy  = hx/hy
332:       hydhx  = hy/hx
333:       temp1  = lambda/(lambda + one)

335: !  Get a pointer to vector data.
336: !    - VecGetArray() returns a pointer to the data array.
337: !    - You MUST call VecRestoreArray() when you no longer need access to
338: !      the array.
339:        call VecGetArray(localX,xx,idx,ierr)

341: !  Get local grid boundaries (for 2-dimensional DA):
342: !    xs, ys   - starting grid indices (no ghost points)
343: !    xm, ym   - widths of local grid (no ghost points)
344: !    gxs, gys - starting grid indices (including ghost points)
345: !    gxm, gym - widths of local grid (including ghost points)

347:        call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,             &
348:      &      PETSC_NULL_INTEGER,ierr)
349:        call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,    &
350:      &      PETSC_NULL_INTEGER,ierr)

352: !  Compute initial guess over the locally owned part of the grid

354:       do 30 j=ys,ys+ym-1
355:         temp = (min(j,my-j-1))*hy
356:         do 40 i=xs,xs+xm-1
357:           row = i - gxs + (j - gys)*gxm + 1
358:           if (i .eq. 0 .or. j .eq. 0 .or. i .eq. mx-1 .or.              &
359:      &        j .eq. my-1) then
360:             xx(idx+row) = 0.d0
361:             continue
362:           endif
363:           xx(idx+row) = temp1*sqrt(min((min(i,mx-i-1))*hx,temp))
364:  40     continue
365:  30   continue

367: !     Restore vector

369:        call VecRestoreArray(localX,xx,idx,ierr)

371: !     Insert values into global vector

373:        call DALocalToGlobal(da,localX,INSERT_VALUES,X,ierr)
374:        return
375:        end

377: ! -------------------------------------------------------------------
378: !
379: !   ComputeFunction - Evaluates nonlinear function, F(x).
380: !
381: !   Input Parameters:
382: !.  X - input vector
383: !
384: !   Output Parameter:
385: !.  F - function vector
386: !
387:       subroutine  ComputeFunction(X,F,ierr)
388:       implicit none

390: !     petsc.h       - base PETSc routines   petscvec.h - vectors
391: !     petscsys.h    - system routines       petscmat.h - matrices
392: !     petscis.h     - index sets            petscksp.h - Krylov subspace methods
393: !     petscviewer.h - viewers               petscpc.h  - preconditioners

395:  #include include/finclude/petsc.h
396:  #include include/finclude/petscis.h
397:  #include include/finclude/petscvec.h
398:  #include include/finclude/petscmat.h
399:  #include include/finclude/petscpc.h
400:  #include include/finclude/petscksp.h
401:  #include include/finclude/petscda.h

403:       Vec              X,F,localX,localF
404:       PetscInt         gys,gxm,gym
405:       PetscOffset      idx,idf
406:       PetscErrorCode ierr
407:       PetscInt i,j,row,mx,my,xs,ys,xm,ym,gxs
408:       double precision two,one,lambda,hx
409:       double precision hy,hxdhy,hydhx,sc
410:       PetscScalar      u,uxx,uyy,xx(1),ff(1)
411:       DA               da
412:       Mat              B
413:       common   /mycommon/ mx,my,B,localX,localF,da

415:       two    = 2.d0
416:       one    = 1.d0
417:       lambda = 6.d0

419:       hx     = one/(mx-1)
420:       hy     = one/(my-1)
421:       sc     = hx*hy*lambda
422:       hxdhy  = hx/hy
423:       hydhx  = hy/hx

425: !  Scatter ghost points to local vector, using the 2-step process
426: !     DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
427: !  By placing code between these two statements, computations can be
428: !  done while messages are in transition.
429: !
430:       call DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)
431:       call DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)

433: !  Get pointers to vector data

435:       call VecGetArray(localX,xx,idx,ierr)
436:       call VecGetArray(localF,ff,idf,ierr)

438: !  Get local grid boundaries

440:       call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,              &
441:      &     PETSC_NULL_INTEGER,ierr)
442:       call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,     &
443:      &     PETSC_NULL_INTEGER,ierr)

445: !  Compute function over the locally owned part of the grid

447:       do 50 j=ys,ys+ym-1

449:         row = (j - gys)*gxm + xs - gxs
450:         do 60 i=xs,xs+xm-1
451:           row = row + 1

453:           if (i .eq. 0 .or. j .eq. 0 .or. i .eq. mx-1 .or.              &
454:      &        j .eq. my-1) then
455:             ff(idf+row) = xx(idx+row)
456:             goto 60
457:           endif
458:           u   = xx(idx+row)
459:           uxx = (two*u - xx(idx+row-1) - xx(idx+row+1))*hydhx
460:           uyy = (two*u - xx(idx+row-gxm) - xx(idx+row+gxm))*hxdhy
461:           ff(idf+row) = uxx + uyy - sc*exp(u)
462:  60     continue
463:  50   continue

465: !  Restore vectors

467:        call VecRestoreArray(localX,xx,idx,ierr)
468:        call VecRestoreArray(localF,ff,idf,ierr)

470: !  Insert values into global vector

472:        call DALocalToGlobal(da,localF,INSERT_VALUES,F,ierr)
473:        return
474:        end

476: ! -------------------------------------------------------------------
477: !
478: !   ComputeJacobian - Evaluates Jacobian matrix.
479: !
480: !   Input Parameters:
481: !   x - input vector
482: !
483: !   Output Parameters:
484: !   jac - Jacobian matrix
485: !   flag - flag indicating matrix structure
486: !
487: !   Notes:
488: !   Due to grid point reordering with DAs, we must always work
489: !   with the local grid points, and then transform them to the new
490: !   global numbering with the 'ltog' mapping (via DAGetGlobalIndices()).
491: !   We cannot work directly with the global numbers for the original
492: !   uniprocessor grid!
493: !
494:       subroutine ComputeJacobian(X,jac,ierr)
495:       implicit none

497: !     petsc.h  - base PETSc routines   petscvec.h - vectors
498: !     petscsys.h    - system routines       petscmat.h - matrices
499: !     petscis.h     - index sets            petscksp.h - Krylov subspace methods
500: !     petscviewer.h - viewers               petscpc.h  - preconditioners

502:  #include include/finclude/petsc.h
503:  #include include/finclude/petscis.h
504:  #include include/finclude/petscvec.h
505:  #include include/finclude/petscmat.h
506:  #include include/finclude/petscpc.h
507:  #include include/finclude/petscksp.h
508:  #include include/finclude/petscda.h

510:       Vec         X
511:       Mat         jac
512:       Vec         localX,localF
513:       DA          da
514:       PetscInt     ltog(1)
515:       PetscOffset idltog,idx
516:       PetscErrorCode ierr
517:       PetscInt nloc,xs,ys,xm,ym
518:       PetscInt gxs,gys,gxm,gym
519:       PetscInt grow(1),i,j
520:       PetscInt row,mx,my,ione
521:       PetscInt col(5),ifive
522:       PetscScalar two,one,lambda
523:       PetscScalar v(5),hx,hy,hxdhy
524:       PetscScalar hydhx,sc,xx(1)
525:       Mat         B
526:       common   /mycommon/ mx,my,B,localX,localF,da

528:       ione   = 1
529:       ifive  = 5
530:       one    = 1.d0
531:       two    = 2.d0
532:       hx     = one/(mx-1)
533:       hy     = one/(my-1)
534:       sc     = hx*hy
535:       hxdhy  = hx/hy
536:       hydhx  = hy/hx
537:       lambda = 6.d0

539: !  Scatter ghost points to local vector, using the 2-step process
540: !     DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
541: !  By placing code between these two statements, computations can be
542: !  done while messages are in transition.

544:       call DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)
545:       call DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)

547: !  Get pointer to vector data

549:       call VecGetArray(localX,xx,idx,ierr)

551: !  Get local grid boundaries

553:       call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,              &
554:      &     PETSC_NULL_INTEGER,ierr)
555:       call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,     &
556:      &                        PETSC_NULL_INTEGER,ierr)

558: !  Get the global node numbers for all local nodes, including ghost points

560:       call DAGetGlobalIndices(da,nloc,ltog,idltog,ierr)

562: !  Compute entries for the locally owned part of the Jacobian.
563: !   - Currently, all PETSc parallel matrix formats are partitioned by
564: !     contiguous chunks of rows across the processors. The 'grow'
565: !     parameter computed below specifies the global row number
566: !     corresponding to each local grid point.
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: !   - Always specify global row and columns of matrix entries.
571: !   - Here, we set all entries for a particular row at once.

573:       do 10 j=ys,ys+ym-1
574:         row = (j - gys)*gxm + xs - gxs
575:         do 20 i=xs,xs+xm-1
576:           row = row + 1
577:           grow(1) = ltog(idltog+row)
578:           if (i .eq. 0 .or. j .eq. 0 .or. i .eq. (mx-1) .or.            &
579:      &        j .eq. (my-1)) then
580:              call MatSetValues(jac,ione,grow,ione,grow,one,             &
581:      &                         INSERT_VALUES,ierr)
582:              go to 20
583:           endif
584:           v(1)   = -hxdhy
585:           col(1) = ltog(idltog+row - gxm)
586:           v(2)   = -hydhx
587:           col(2) = ltog(idltog+row - 1)
588:           v(3)   = two*(hydhx + hxdhy) - sc*lambda*exp(xx(idx+row))
589:           col(3) = grow(1)
590:           v(4)   = -hydhx
591:           col(4) = ltog(idltog+row + 1)
592:           v(5)   = -hxdhy
593:           col(5) = ltog(idltog+row + gxm)
594:           call MatSetValues(jac,ione,grow,ifive,col,v,INSERT_VALUES,       &
595:      &                      ierr)
596:  20     continue
597:  10   continue

599: !  Assemble matrix, using the 2-step process:
600: !    MatAssemblyBegin(), MatAssemblyEnd().
601: !  By placing code between these two statements, computations can be
602: !  done while messages are in transition.

604:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
605:       call VecRestoreArray(localX,xx,idx,ierr)
606:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
607:       return
608:       end


611: ! -------------------------------------------------------------------
612: !
613: !   MyMult - user provided matrix multiply
614: !
615: !   Input Parameters:
616: !.  X - input vector
617: !
618: !   Output Parameter:
619: !.  F - function vector
620: !
621:       subroutine  MyMult(J,X,F,ierr)
622:       implicit none
623:       Mat     J,B
624:       Vec     X,F
625:       PetscErrorCode ierr
626:       PetscInt mx,my
627:       DA      da
628:       Vec     localX,localF

630:       common   /mycommon/ mx,my,B,localX,localF,da
631: !
632: !       Here we use the actual formed matrix B; users would
633: !     instead write their own matrix vector product routine
634: !
635:       call MatMult(B,X,F,ierr)
636:       return
637:       end