Actual source code: ex1f.F

  1: ! "$Id: ex1f.F,v 1.56 2001/08/10 16:48:34 balay Exp $";
  2: !
  3: !  Description: This example solves a nonlinear system on 1 processor with SNES.
  4: !  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
  5: !  domain.  The command line options include:
  6: !    -par <parameter>, where <parameter> indicates the nonlinearity of the problem
  7: !       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
  8: !    -mx <xg>, where <xg> = number of grid points in the x-direction
  9: !    -my <yg>, where <yg> = number of grid points in the y-direction
 10: !
 11: !/*T
 12: !  Concepts: SNES^sequential Bratu example
 13: !  Processors: 1
 14: !T*/
 15: !
 16: !  --------------------------------------------------------------------------
 17: !
 18: !  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 19: !  the partial differential equation
 20: !
 21: !          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 22: !
 23: !  with boundary conditions
 24: !
 25: !           u = 0  for  x = 0, x = 1, y = 0, y = 1.
 26: !
 27: !  A finite difference approximation with the usual 5-point stencil
 28: !  is used to discretize the boundary value problem to obtain a nonlinear
 29: !  system of equations.
 30: !
 31: !  The parallel version of this code is snes/examples/tutorials/ex5f.F
 32: !
 33: !  --------------------------------------------------------------------------

 35:       program main
 36:       implicit none

 38: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 39: !                    Include files
 40: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41: !
 42: !  The following include statements are generally used in SNES Fortran
 43: !  programs:
 44: !     petsc.h       - base PETSc routines
 45: !     petscvec.h    - vectors
 46: !     petscmat.h    - matrices
 47: !     petscksp.h    - Krylov subspace methods
 48: !     petscpc.h     - preconditioners
 49: !     petscsnes.h   - SNES interface
 50: !  In addition, we need the following for use of PETSc drawing routines
 51: !     petscdraw.h   - drawing routines
 52: !  Other include statements may be needed if using additional PETSc
 53: !  routines in a Fortran program, e.g.,
 54: !     petscviewer.h - viewers
 55: !     petscis.h     - index sets
 56: !
 57:  #include include/finclude/petsc.h
 58:  #include include/finclude/petscvec.h
 59:  #include include/finclude/petscis.h
 60:  #include include/finclude/petscdraw.h
 61:  #include include/finclude/petscmat.h
 62:  #include include/finclude/petscksp.h
 63:  #include include/finclude/petscpc.h
 64:  #include include/finclude/petscsnes.h
 65: !
 66: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67: !                   Variable declarations
 68: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69: !
 70: !  Variables:
 71: !     snes        - nonlinear solver
 72: !     x, r        - solution, residual vectors
 73: !     J           - Jacobian matrix
 74: !     its         - iterations for convergence
 75: !     matrix_free - flag - 1 indicates matrix-free version
 76: !     lambda      - nonlinearity parameter
 77: !     draw        - drawing context
 78: !
 79:       SNES               snes
 80:       Vec                x,r
 81:       PetscDraw               draw
 82:       Mat                J
 83:       integer            its,matrix_free,flg,N,ierr
 84:       integer            mx,my,size,rank
 85:       integer            fd_coloring
 86:       PetscReal          lambda_max,lambda_min,lambda
 87:       MatFDColoring      fdcoloring
 88:       ISColoring         iscoloring
 89:       MatStructure       str
 90: 
 91:       PetscScalar        lx_v(0:1)
 92:       PetscOffset        lx_i

 94: !  Store parameters in common block

 96:       common /params/ lambda,mx,my

 98: !  Note: Any user-defined Fortran routines (such as FormJacobian)
 99: !  MUST be declared as external.

101:       external FormFunction,FormInitialGuess,FormJacobian

103: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: !  Initialize program
105: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

107:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
108:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
109:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

111:       if (size .ne. 1) then
112:          if (rank .eq. 0) then
113:             write(6,*) 'This is a uniprocessor example only!'
114:          endif
115:          SETERRQ(1,' ',ierr)
116:       endif

118: !  Initialize problem parameters

120:       lambda_max = 6.81
121:       lambda_min = 0.0
122:       lambda     = 6.0
123:       mx         = 4
124:       my         = 4
125:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr)
126:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-my',my,flg,ierr)
127:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',lambda,      &
128:      &     flg,ierr)
129:       if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
130:          if (rank .eq. 0) write(6,*) 'Lambda is out of range'
131:          SETERRQ(1,' ',ierr)
132:       endif
133:       N       = mx*my

135: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: !  Create nonlinear solver context
137: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

139:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

141: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: !  Create vector data structures; set function evaluation routine
143: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

145:       call VecCreate(PETSC_COMM_WORLD,x,ierr)
146:       call VecSetSizes(x,PETSC_DECIDE,N,ierr)
147:       call VecSetFromOptions(x,ierr)
148:       call VecDuplicate(x,r,ierr)

150: !  Set function evaluation routine and vector.  Whenever the nonlinear
151: !  solver needs to evaluate the nonlinear function, it will call this
152: !  routine.
153: !   - Note that the final routine argument is the user-defined
154: !     context that provides application-specific data for the
155: !     function evaluation routine.

157:       call SNESSetFunction(snes,r,FormFunction,PETSC_NULL_OBJECT,ierr)

159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: !  Create matrix data structure; set Jacobian evaluation routine
161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

163: !  Create matrix. Here we only approximately preallocate storage space
164: !  for the Jacobian.  See the users manual for a discussion of better
165: !  techniques for preallocating matrix memory.

167:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_mf',         &
168:      &     matrix_free,ierr)
169:       if (matrix_free .eq. 0) then
170:         call MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,PETSC_NULL_INTEGER, &
171:      &        J,ierr)
172:       endif

174: !
175: !     This option will cause the Jacobian to be computed via finite differences
176: !    efficiently using a coloring of the columns of the matrix.
177: !
178:       fd_coloring = 0
179:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_fd_coloring',     &
180:      &                    fd_coloring,ierr)
181:       if (fd_coloring .eq. 1) then
182: 
183: !
184: !      This initializes the nonzero structure of the Jacobian. This is artificial
185: !      because clearly if we had a routine to compute the Jacobian we won't need
186: !      to use finite differences.
187: !
188:         call FormJacobian(snes,x,J,J,str,0,ierr)
189: !
190: !       Color the matrix, i.e. determine groups of columns that share no common
191: !      rows. These columns in the Jacobian can all be computed simulataneously.
192: !
193:         call MatGetColoring(J,MATCOLORING_NATURAL,iscoloring,ierr)
194: 
195: !
196: !       Create the data structure that SNESDefaultComputeJacobianColor() uses
197: !       to compute the actual Jacobians via finite differences.
198: !
199:         call MatFDColoringCreate(J,iscoloring,fdcoloring,ierr)
200:         call MatFDColoringSetFunctionSNES(fdcoloring,FormFunction,        &
201:      &                                PETSC_NULL_OBJECT,ierr)
202:         call MatFDColoringSetFromOptions(fdcoloring,ierr)
203: !
204: !        Tell SNES to use the routine SNESDefaultComputeJacobianColor()
205: !      to compute Jacobians.
206: !
207:         call SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,    &
208:      &                     fdcoloring,ierr)
209:         call ISColoringDestroy(iscoloring,ierr)

211:       else if (matrix_free .eq. 0) then

213: !  Set Jacobian matrix data structure and default Jacobian evaluation
214: !  routine.  Whenever the nonlinear solver needs to compute the
215: !  Jacobian matrix, it will call this routine.
216: !   - Note that the final routine argument is the user-defined
217: !     context that provides application-specific data for the
218: !     Jacobian evaluation routine.
219: !   - The user can override with:
220: !      -snes_fd : default finite differencing approximation of Jacobian
221: !      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
222: !                 (unless user explicitly sets preconditioner)
223: !      -snes_mf_operator : form preconditioning matrix as set by the user,
224: !                          but use matrix-free approx for Jacobian-vector
225: !                          products within Newton-Krylov method
226: !
227:         call SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL_OBJECT,   &
228:      &        ierr)
229:       endif

231: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232: !  Customize nonlinear solver; set runtime options
233: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

235: !  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)

237:       call SNESSetFromOptions(snes,ierr)

239: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
240: !  Evaluate initial guess; then solve nonlinear system.
241: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

243: !  Note: The user should initialize the vector, x, with the initial guess
244: !  for the nonlinear solver prior to calling SNESSolve().  In particular,
245: !  to employ an initial guess of zero, the user should explicitly set
246: !  this vector to zero by calling VecSet().

248:       call FormInitialGuess(x,ierr)
249:       call SNESSolve(snes,x,ierr)
250:       call SNESGetIterationNumber(snes,its,ierr);
251:       if (rank .eq. 0) then
252:          write(6,100) its
253:       endif
254:   100 format('Number of Newton iterations = ',i1)

256: !  PetscDraw contour plot of solution

258: !      call PetscDrawOpenX(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,'Solution',
259: !                   300,0,300,300,draw,ierr)
260:       call PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,          &
261:      &     'Solution',300,0,300,300,draw,ierr)
262:       call PetscDrawSetType(draw,PETSC_DRAW_X,ierr)

264:       call VecGetArray(x,lx_v,lx_i,ierr)
265:       call PetscDrawTensorContour(draw,mx,my,PETSC_NULL_DOUBLE,              &
266:      &     PETSC_NULL_DOUBLE,lx_v(lx_i+1),ierr)
267:       call VecRestoreArray(x,lx_v,lx_i,ierr)

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

274:       if (matrix_free .eq. 0) call MatDestroy(J,ierr)
275:       if (fd_coloring .ne. 0) call MatFDColoringDestroy(fdcoloring,ierr)

277:       call VecDestroy(x,ierr)
278:       call VecDestroy(r,ierr)
279:       call SNESDestroy(snes,ierr)
280:       call PetscDrawDestroy(draw,ierr)
281:       call PetscFinalize(ierr)
282:       end

284: ! ---------------------------------------------------------------------
285: !
286: !  FormInitialGuess - Forms initial approximation.
287: !
288: !  Input Parameter:
289: !  X - vector
290: !
291: !  Output Parameters:
292: !  X - vector
293: !  ierr - error code
294: !
295: !  Notes:
296: !  This routine serves as a wrapper for the lower-level routine
297: !  "ApplicationInitialGuess", where the actual computations are
298: !  done using the standard Fortran style of treating the local
299: !  vector data as a multidimensional array over the local mesh.
300: !  This routine merely accesses the local vector data via
301: !  VecGetArray() and VecRestoreArray().
302: !
303:       subroutine FormInitialGuess(X,ierr)
304:       implicit none

306:  #include include/finclude/petsc.h
307:  #include include/finclude/petscvec.h
308:  #include include/finclude/petscmat.h
309:  #include include/finclude/petscsnes.h

311: !  Input/output variables:
312:       Vec           X
313:       integer       ierr

315: !  Declarations for use with local arrays:
316:       PetscScalar   lx_v(0:1)
317:       PetscOffset   lx_i

319:       0

321: !  Get a pointer to vector data.
322: !    - For default PETSc vectors, VecGetArray() returns a pointer to
323: !      the data array.  Otherwise, the routine is implementation dependent.
324: !    - You MUST call VecRestoreArray() when you no longer need access to
325: !      the array.
326: !    - Note that the Fortran interface to VecGetArray() differs from the
327: !      C version.  See the users manual for details.

329:       call VecGetArray(X,lx_v,lx_i,ierr)

331: !  Compute initial guess

333:       call ApplicationInitialGuess(lx_v(lx_i),ierr)

335: !  Restore vector

337:       call VecRestoreArray(X,lx_v,lx_i,ierr)

339:       return
340:       end

342: ! ---------------------------------------------------------------------
343: !
344: !  ApplicationInitialGuess - Computes initial approximation, called by
345: !  the higher level routine FormInitialGuess().
346: !
347: !  Input Parameter:
348: !  x - local vector data
349: !
350: !  Output Parameters:
351: !  f - local vector data, f(x)
352: !  ierr - error code
353: !
354: !  Notes:
355: !  This routine uses standard Fortran-style computations over a 2-dim array.
356: !
357:       subroutine ApplicationInitialGuess(x,ierr)

359:       implicit none

361: !  Common blocks:
362:       PetscReal   lambda
363:       integer     mx,my
364:       common      /params/ lambda,mx,my

366: !  Input/output variables:
367:       PetscScalar x(mx,my)
368:       integer     ierr

370: !  Local variables:
371:       integer     i,j,hxdhy,hydhx
372:       PetscScalar temp1,temp,hx,hy,sc,one

374: !  Set parameters

376:       0
377:       one    = 1.0
378:       hx     = one/(dble(mx-1))
379:       hy     = one/(dble(my-1))
380:       sc     = hx*hy*lambda
381:       hxdhy  = hx/hy
382:       hydhx  = hy/hx
383:       temp1  = lambda/(lambda + one)

385:       do 20 j=1,my
386:          temp = dble(min(j-1,my-j))*hy
387:          do 10 i=1,mx
388:             if (i .eq. 1 .or. j .eq. 1                                  &
389:      &             .or. i .eq. mx .or. j .eq. my) then
390:               x(i,j) = 0.0
391:             else
392:               x(i,j) = temp1 *                                          &
393:      &          sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
394:             endif
395:  10      continue
396:  20   continue

398:       return
399:       end

401: ! ---------------------------------------------------------------------
402: !
403: !  FormFunction - Evaluates nonlinear function, F(x).
404: !
405: !  Input Parameters:
406: !  snes  - the SNES context
407: !  X     - input vector
408: !  dummy - optional user-defined context, as set by SNESSetFunction()
409: !          (not used here)
410: !
411: !  Output Parameter:
412: !  F     - vector with newly computed function
413: !
414: !  Notes:
415: !  This routine serves as a wrapper for the lower-level routine
416: !  "ApplicationFunction", where the actual computations are
417: !  done using the standard Fortran style of treating the local
418: !  vector data as a multidimensional array over the local mesh.
419: !  This routine merely accesses the local vector data via
420: !  VecGetArray() and VecRestoreArray().
421: !
422:       subroutine FormFunction(snes,X,F,dummy,ierr)
423:       implicit none

425:  #include include/finclude/petsc.h
426:  #include include/finclude/petscvec.h
427:  #include include/finclude/petscsnes.h

429: !  Input/output variables:
430:       SNES              snes
431:       Vec               X,F
432:       PetscFortranAddr  dummy
433:       integer           ierr

435: !  Common blocks:
436:       PetscReal         lambda
437:       integer           mx,my
438:       common            /params/ lambda,mx,my

440: !  Declarations for use with local arrays:
441:       PetscScalar       lx_v(0:1),lf_v(0:1)
442:       PetscOffset       lx_i,lf_i

444: !  Get pointers to vector data.
445: !    - For default PETSc vectors, VecGetArray() returns a pointer to
446: !      the data array.  Otherwise, the routine is implementation dependent.
447: !    - You MUST call VecRestoreArray() when you no longer need access to
448: !      the array.
449: !    - Note that the Fortran interface to VecGetArray() differs from the
450: !      C version.  See the Fortran chapter of the users manual for details.

452:       call VecGetArray(X,lx_v,lx_i,ierr)
453:       call VecGetArray(F,lf_v,lf_i,ierr)

455: !  Compute function

457:       call ApplicationFunction(lx_v(lx_i),lf_v(lf_i),ierr)

459: !  Restore vectors

461:       call VecRestoreArray(X,lx_v,lx_i,ierr)
462:       call VecRestoreArray(F,lf_v,lf_i,ierr)

464:       call PetscLogFlops(11*mx*my,ierr)

466:       return
467:       end

469: ! ---------------------------------------------------------------------
470: !
471: !  ApplicationFunction - Computes nonlinear function, called by
472: !  the higher level routine FormFunction().
473: !
474: !  Input Parameter:
475: !  x    - local vector data
476: !
477: !  Output Parameters:
478: !  f    - local vector data, f(x)
479: !  ierr - error code
480: !
481: !  Notes:
482: !  This routine uses standard Fortran-style computations over a 2-dim array.
483: !
484:       subroutine ApplicationFunction(x,f,ierr)

486:       implicit none

488: !  Common blocks:
489:       PetscReal      lambda
490:       integer        mx,my
491:       common         /params/ lambda,mx,my

493: !  Input/output variables:
494:       PetscScalar    x(mx,my),f(mx,my)
495:       integer        ierr

497: !  Local variables:
498:       PetscScalar    two,one,hx,hy,hxdhy,hydhx,sc
499:       PetscScalar    u,uxx,uyy
500:       integer        i,j

502:       0
503:       one    = 1.0
504:       two    = 2.0
505:       hx     = one/dble(mx-1)
506:       hy     = one/dble(my-1)
507:       sc     = hx*hy*lambda
508:       hxdhy  = hx/hy
509:       hydhx  = hy/hx

511: !  Compute function

513:       do 20 j=1,my
514:          do 10 i=1,mx
515:             if (i .eq. 1 .or. j .eq. 1                                  &
516:      &             .or. i .eq. mx .or. j .eq. my) then
517:                f(i,j) = x(i,j)
518:             else
519:                u = x(i,j)
520:                uxx = hydhx * (two*u                                     &
521:      &                - x(i-1,j) - x(i+1,j))
522:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
523:                f(i,j) = uxx + uyy - sc*exp(u)
524:             endif
525:  10      continue
526:  20   continue

528:       return
529:       end

531: ! ---------------------------------------------------------------------
532: !
533: !  FormJacobian - Evaluates Jacobian matrix.
534: !
535: !  Input Parameters:
536: !  snes    - the SNES context
537: !  x       - input vector
538: !  dummy   - optional user-defined context, as set by SNESSetJacobian()
539: !            (not used here)
540: !
541: !  Output Parameters:
542: !  jac      - Jacobian matrix
543: !  jac_prec - optionally different preconditioning matrix (not used here)
544: !  flag     - flag indicating matrix structure
545: !
546: !  Notes:
547: !  This routine serves as a wrapper for the lower-level routine
548: !  "ApplicationJacobian", where the actual computations are
549: !  done using the standard Fortran style of treating the local
550: !  vector data as a multidimensional array over the local mesh.
551: !  This routine merely accesses the local vector data via
552: !  VecGetArray() and VecRestoreArray().
553: !
554:       subroutine FormJacobian(snes,X,jac,jac_prec,flag,dummy,ierr)
555:       implicit none

557:  #include include/finclude/petsc.h
558:  #include include/finclude/petscvec.h
559:  #include include/finclude/petscmat.h
560:  #include include/finclude/petscpc.h
561:  #include include/finclude/petscsnes.h

563: !  Input/output variables:
564:       SNES          snes
565:       Vec           X
566:       Mat           jac,jac_prec
567:       MatStructure  flag
568:       integer       dummy,ierr

570: !  Common blocks:
571:       PetscReal     lambda
572:       integer       mx,my
573:       common        /params/ lambda,mx,my

575: !  Declarations for use with local array:
576:       PetscScalar   lx_v(0:1)
577:       PetscOffset   lx_i

579: !  Get a pointer to vector data

581:       call VecGetArray(X,lx_v,lx_i,ierr)

583: !  Compute Jacobian entries

585:       call ApplicationJacobian(lx_v(lx_i),jac,jac_prec,ierr)

587: !  Restore vector

589:       call VecRestoreArray(X,lx_v,lx_i,ierr)

591: !  Assemble matrix

593:       call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
594:       call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)

596: !  Set flag to indicate that the Jacobian matrix retains an identical
597: !  nonzero structure throughout all nonlinear iterations (although the
598: !  values of the entries change). Thus, we can save some work in setting
599: !  up the preconditioner (e.g., no need to redo symbolic factorization for
600: !  ILU/ICC preconditioners).
601: !   - If the nonzero structure of the matrix is different during
602: !     successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
603: !     must be used instead.  If you are unsure whether the matrix
604: !     structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
605: !   - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
606: !     believes your assertion and does not check the structure
607: !     of the matrix.  If you erroneously claim that the structure
608: !     is the same when it actually is not, the new preconditioner
609: !     will not function correctly.  Thus, use this optimization
610: !     feature with caution!

612:       flag = SAME_NONZERO_PATTERN

614:       return
615:       end

617: ! ---------------------------------------------------------------------
618: !
619: !  ApplicationJacobian - Computes Jacobian matrix, called by
620: !  the higher level routine FormJacobian().
621: !
622: !  Input Parameters:
623: !  x        - local vector data
624: !
625: !  Output Parameters:
626: !  jac      - Jacobian matrix
627: !  jac_prec - optionally different preconditioning matrix (not used here)
628: !  ierr     - error code
629: !
630: !  Notes:
631: !  This routine uses standard Fortran-style computations over a 2-dim array.
632: !
633:       subroutine ApplicationJacobian(x,jac,jac_prec,ierr)
634:       implicit none

636:  #include include/finclude/petsc.h
637:  #include include/finclude/petscvec.h
638:  #include include/finclude/petscmat.h
639:  #include include/finclude/petscpc.h
640:  #include include/finclude/petscsnes.h

642: !  Common blocks:
643:       PetscReal    lambda
644:       integer      mx,my
645:       common       /params/ lambda,mx,my

647: !  Input/output variables:
648:       PetscScalar  x(mx,my)
649:       Mat          jac,jac_prec
650:       integer      ierr

652: !  Local variables:
653:       integer      i,j,row,col(5)
654:       PetscScalar  two,one, hx,hy,hxdhy,hydhx,sc,v(5)

656: !  Set parameters

658:       one    = 1.0
659:       two    = 2.0
660:       hx     = one/dble(mx-1)
661:       hy     = one/dble(my-1)
662:       sc     = hx*hy
663:       hxdhy  = hx/hy
664:       hydhx  = hy/hx

666: !  Compute entries of the Jacobian matrix
667: !   - Here, we set all entries for a particular row at once.
668: !   - Note that MatSetValues() uses 0-based row and column numbers
669: !     in Fortran as well as in C.

671:       do 20 j=1,my
672:          row = (j-1)*mx - 1
673:          do 10 i=1,mx
674:             row = row + 1
675: !           boundary points
676:             if (i .eq. 1 .or. j .eq. 1                                  &
677:      &             .or. i .eq. mx .or. j .eq. my) then
678:                call MatSetValues(jac_prec,1,row,1,row,one,              &
679:      &                           INSERT_VALUES,ierr)
680: !           interior grid points
681:             else
682:                v(1) = -hxdhy
683:                v(2) = -hydhx
684:                v(3) = two*(hydhx + hxdhy)                               &
685:      &                  - sc*lambda*exp(x(i,j))
686:                v(4) = -hydhx
687:                v(5) = -hxdhy
688:                col(1) = row - mx
689:                col(2) = row - 1
690:                col(3) = row
691:                col(4) = row + 1
692:                col(5) = row + mx
693:                call MatSetValues(jac_prec,1,row,5,col,v,                &
694:      &                           INSERT_VALUES,ierr)
695:             endif
696:  10      continue
697:  20   continue

699:       return
700:       end