Actual source code: ex1f.F

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

 34:       program main
 35:       implicit none

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

 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(PETSC_COMM_SELF,1,' ',ierr)
116:       endif

118: !  Initialize problem parameters
119:       i5 = 5
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_OBJECT,PETSC_NULL_CHARACTER,   &
126:      &                        '-mx',mx,flg,ierr)
127:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,   &
128:      &                        '-my',my,flg,ierr)
129:       call PetscOptionsGetReal(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,  &
130:      &                         '-par',lambda,flg,ierr)
131:       if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
132:          if (rank .eq. 0) write(6,*) 'Lambda is out of range'
133:          SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
134:       endif
135:       N       = mx*my

137: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: !  Create nonlinear solver context
139: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

141:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

143: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: !  Create vector data structures; set function evaluation routine
145: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

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

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

161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162: !  Create matrix data structure; set Jacobian evaluation routine
163: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

169:       call PetscOptionsHasName(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,                    &
170:      &                         '-snes_mf',matrix_free,ierr)
171:       if (.not. matrix_free) then
172:         call MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,i5,PETSC_NULL_INTEGER, &
173:      &        J,ierr)
174:       endif

176: !
177: !     This option will cause the Jacobian to be computed via finite differences
178: !    efficiently using a coloring of the columns of the matrix.
179: !
180:       fd_coloring = .false.
181:       call PetscOptionsHasName(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,       &
182:      &                         '-snes_fd_coloring',fd_coloring,ierr)
183:       if (fd_coloring) then

185: !
186: !      This initializes the nonzero structure of the Jacobian. This is artificial
187: !      because clearly if we had a routine to compute the Jacobian we won't need
188: !      to use finite differences.
189: !
190:         call FormJacobian(snes,x,J,J,0,ierr)
191: !
192: !       Color the matrix, i.e. determine groups of columns that share no common
193: !      rows. These columns in the Jacobian can all be computed simulataneously.
194: !
195:         call MatColoringCreate(J,mc,ierr)
196:         call MatColoringSetType(mc,MATCOLORINGNATURAL,ierr)
197:         call MatColoringSetFromOptions(mc,ierr)
198:         call MatColoringApply(mc,iscoloring,ierr)
199:         call MatColoringDestroy(mc,ierr)
200: !
201: !       Create the data structure that SNESComputeJacobianDefaultColor() uses
202: !       to compute the actual Jacobians via finite differences.
203: !
204:         call MatFDColoringCreate(J,iscoloring,fdcoloring,ierr)
205:         call MatFDColoringSetFunction(fdcoloring,FormFunction,                &
206:      &                                PETSC_NULL_OBJECT,ierr)
207:         call MatFDColoringSetFromOptions(fdcoloring,ierr)
208:         call MatFDColoringSetUp(J,iscoloring,fdcoloring,ierr)
209: !
210: !        Tell SNES to use the routine SNESComputeJacobianDefaultColor()
211: !      to compute Jacobians.
212: !
213:         call SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,    &
214:      &                     fdcoloring,ierr)
215:         call ISColoringDestroy(iscoloring,ierr)

217:       else if (.not. matrix_free) then

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

237: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238: !  Customize nonlinear solver; set runtime options
239: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

243:       call SNESSetFromOptions(snes,ierr)

245: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
246: !  Evaluate initial guess; then solve nonlinear system.
247: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

254:       call FormInitialGuess(x,ierr)
255:       call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
256:       call SNESGetIterationNumber(snes,its,ierr);
257:       if (rank .eq. 0) then
258:          write(6,100) its
259:       endif
260:   100 format('Number of SNES iterations = ',i1)

262: !  PetscDraw contour plot of solution

264:       call PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,          &
265:      &     'Solution',300,0,300,300,draw,ierr)
266:       call PetscDrawSetFromOptions(draw,ierr)

268:       call VecGetArrayRead(x,lx_v,lx_i,ierr)
269:       call PetscDrawTensorContour(draw,mx,my,PETSC_NULL_REAL,              &
270:      &     PETSC_NULL_REAL,lx_v(lx_i+1),ierr)
271:       call VecRestoreArrayRead(x,lx_v,lx_i,ierr)

273: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
274: !  Free work space.  All PETSc objects should be destroyed when they
275: !  are no longer needed.
276: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

278:       if (.not. matrix_free) call MatDestroy(J,ierr)
279:       if (fd_coloring) call MatFDColoringDestroy(fdcoloring,ierr)

281:       call VecDestroy(x,ierr)
282:       call VecDestroy(r,ierr)
283:       call SNESDestroy(snes,ierr)
284:       call PetscDrawDestroy(draw,ierr)
285:       call PetscFinalize(ierr)
286:       end

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

310: #include <petsc/finclude/petscsys.h>
311: #include <petsc/finclude/petscvec.h>
312: #include <petsc/finclude/petscmat.h>
313: #include <petsc/finclude/petscsnes.h>

315: !  Input/output variables:
316:       Vec           X
317:       PetscErrorCode    ierr

319: !  Declarations for use with local arrays:
320:       PetscScalar   lx_v(0:1)
321:       PetscOffset   lx_i

323:       0

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

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

335: !  Compute initial guess

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

339: !  Restore vector

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

343:       return
344:       end

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

363:       implicit none

365: !  Common blocks:
366:       PetscReal   lambda
367:       PetscInt     mx,my
368:       common      /params/ lambda,mx,my

370: !  Input/output variables:
371:       PetscScalar x(mx,my)
372:       PetscErrorCode     ierr

374: !  Local variables:
375:       PetscInt     i,j
376:       PetscReal temp1,temp,hx,hy,one

378: !  Set parameters

380:       0
381:       one    = 1.0
382:       hx     = one/(mx-1)
383:       hy     = one/(my-1)
384:       temp1  = lambda/(lambda + one)

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

399:       return
400:       end

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

426: #include <petsc/finclude/petscsys.h>
427: #include <petsc/finclude/petscvec.h>
428: #include <petsc/finclude/petscsnes.h>

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

436: !  Common blocks:
437:       PetscReal         lambda
438:       PetscInt          mx,my
439:       common            /params/ lambda,mx,my

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

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

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

456: !  Compute function

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

460: !  Restore vectors

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

465:       call PetscLogFlops(11.0d0*mx*my,ierr)

467:       return
468:       end

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

487:       implicit none

489: !  Common blocks:
490:       PetscReal      lambda
491:       PetscInt        mx,my
492:       common         /params/ lambda,mx,my

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

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

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

513: !  Compute function

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

530:       return
531:       end

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

559: #include <petsc/finclude/petscsys.h>
560: #include <petsc/finclude/petscvec.h>
561: #include <petsc/finclude/petscmat.h>
562: #include <petsc/finclude/petscpc.h>
563: #include <petsc/finclude/petscsnes.h>

565: !  Input/output variables:
566:       SNES          snes
567:       Vec           X
568:       Mat           jac,jac_prec
569:       PetscErrorCode      ierr
570:       integer dummy

572: !  Common blocks:
573:       PetscReal     lambda
574:       PetscInt       mx,my
575:       common        /params/ lambda,mx,my

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

581: !  Get a pointer to vector data

583:       call VecGetArrayRead(X,lx_v,lx_i,ierr)

585: !  Compute Jacobian entries

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

589: !  Restore vector

591:       call VecRestoreArrayRead(X,lx_v,lx_i,ierr)

593: !  Assemble matrix

595:       call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
596:       call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)


599:       return
600:       end

602: ! ---------------------------------------------------------------------
603: !
604: !  ApplicationJacobian - Computes Jacobian matrix, called by
605: !  the higher level routine FormJacobian().
606: !
607: !  Input Parameters:
608: !  x        - local vector data
609: !
610: !  Output Parameters:
611: !  jac      - Jacobian matrix
612: !  jac_prec - optionally different preconditioning matrix (not used here)
613: !  ierr     - error code
614: !
615: !  Notes:
616: !  This routine uses standard Fortran-style computations over a 2-dim array.
617: !
618:       subroutine ApplicationJacobian(x,jac,jac_prec,ierr)
619:       implicit none

621: #include <petsc/finclude/petscsys.h>
622: #include <petsc/finclude/petscvec.h>
623: #include <petsc/finclude/petscmat.h>
624: #include <petsc/finclude/petscpc.h>
625: #include <petsc/finclude/petscsnes.h>

627: !  Common blocks:
628:       PetscReal    lambda
629:       PetscInt      mx,my
630:       common       /params/ lambda,mx,my

632: !  Input/output variables:
633:       PetscScalar  x(mx,my)
634:       Mat          jac,jac_prec
635:       PetscErrorCode      ierr

637: !  Local variables:
638:       PetscInt      i,j,row(1),col(5),i1,i5
639:       PetscScalar  two,one, hx,hy
640:       PetscScalar  hxdhy,hydhx,sc,v(5)

642: !  Set parameters

644:       i1     = 1
645:       i5     = 5
646:       one    = 1.0
647:       two    = 2.0
648:       hx     = one/(mx-1)
649:       hy     = one/(my-1)
650:       sc     = hx*hy
651:       hxdhy  = hx/hy
652:       hydhx  = hy/hx

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

659:       do 20 j=1,my
660:          row(1) = (j-1)*mx - 1
661:          do 10 i=1,mx
662:             row(1) = row(1) + 1
663: !           boundary points
664:             if (i .eq. 1 .or. j .eq. 1                                  &
665:      &             .or. i .eq. mx .or. j .eq. my) then
666:                call MatSetValues(jac_prec,i1,row,i1,row,one,              &
667:      &                           INSERT_VALUES,ierr)
668: !           interior grid points
669:             else
670:                v(1) = -hxdhy
671:                v(2) = -hydhx
672:                v(3) = two*(hydhx + hxdhy)                               &
673:      &                  - sc*lambda*exp(x(i,j))
674:                v(4) = -hydhx
675:                v(5) = -hxdhy
676:                col(1) = row(1) - mx
677:                col(2) = row(1) - 1
678:                col(3) = row(1)
679:                col(4) = row(1) + 1
680:                col(5) = row(1) + mx
681:                call MatSetValues(jac_prec,i1,row,i5,col,v,                &
682:      &                           INSERT_VALUES,ierr)
683:             endif
684:  10      continue
685:  20   continue

687:       return
688:       end