Actual source code: ex1f.F
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: ! petsc.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 include/finclude/petsc.h
57: #include include/finclude/petscvec.h
58: #include include/finclude/petscis.h
59: #include include/finclude/petscdraw.h
60: #include include/finclude/petscmat.h
61: #include include/finclude/petscksp.h
62: #include include/finclude/petscpc.h
63: #include include/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: Vec x,r
80: PetscDraw draw
81: Mat J
82: PetscTruth matrix_free,flg,fd_coloring
83: PetscErrorCode ierr
84: PetscInt its,N, mx,my,i5
85: PetscMPIInt size,rank
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
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_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,i5,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,PETSC_NULL_OBJECT,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: PetscErrorCode 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: PetscInt mx,my
364: common /params/ lambda,mx,my
366: ! Input/output variables:
367: PetscScalar x(mx,my)
368: PetscErrorCode ierr
370: ! Local variables:
371: PetscInt 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: PetscErrorCode ierr
435: ! Common blocks:
436: PetscReal lambda
437: PetscInt 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: PetscInt mx,my
491: common /params/ lambda,mx,my
493: ! Input/output variables:
494: PetscScalar x(mx,my),f(mx,my)
495: PetscErrorCode ierr
497: ! Local variables:
498: PetscScalar two,one,hx,hy
499: PetscScalar hxdhy,hydhx,sc
500: PetscScalar u,uxx,uyy
501: PetscInt i,j
503: 0
504: one = 1.0
505: two = 2.0
506: hx = one/dble(mx-1)
507: hy = one/dble(my-1)
508: sc = hx*hy*lambda
509: hxdhy = hx/hy
510: hydhx = hy/hx
512: ! Compute function
514: do 20 j=1,my
515: do 10 i=1,mx
516: if (i .eq. 1 .or. j .eq. 1 &
517: & .or. i .eq. mx .or. j .eq. my) then
518: f(i,j) = x(i,j)
519: else
520: u = x(i,j)
521: uxx = hydhx * (two*u &
522: & - x(i-1,j) - x(i+1,j))
523: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
524: f(i,j) = uxx + uyy - sc*exp(u)
525: endif
526: 10 continue
527: 20 continue
529: return
530: end
532: ! ---------------------------------------------------------------------
533: !
534: ! FormJacobian - Evaluates Jacobian matrix.
535: !
536: ! Input Parameters:
537: ! snes - the SNES context
538: ! x - input vector
539: ! dummy - optional user-defined context, as set by SNESSetJacobian()
540: ! (not used here)
541: !
542: ! Output Parameters:
543: ! jac - Jacobian matrix
544: ! jac_prec - optionally different preconditioning matrix (not used here)
545: ! flag - flag indicating matrix structure
546: !
547: ! Notes:
548: ! This routine serves as a wrapper for the lower-level routine
549: ! "ApplicationJacobian", where the actual computations are
550: ! done using the standard Fortran style of treating the local
551: ! vector data as a multidimensional array over the local mesh.
552: ! This routine merely accesses the local vector data via
553: ! VecGetArray() and VecRestoreArray().
554: !
555: subroutine FormJacobian(snes,X,jac,jac_prec,flag,dummy,ierr)
556: implicit none
558: #include include/finclude/petsc.h
559: #include include/finclude/petscvec.h
560: #include include/finclude/petscmat.h
561: #include include/finclude/petscpc.h
562: #include include/finclude/petscsnes.h
564: ! Input/output variables:
565: SNES snes
566: Vec X
567: Mat jac,jac_prec
568: MatStructure flag
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 VecGetArray(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 VecRestoreArray(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)
598: ! Set flag to indicate that the Jacobian matrix retains an identical
599: ! nonzero structure throughout all nonlinear iterations (although the
600: ! values of the entries change). Thus, we can save some work in setting
601: ! up the preconditioner (e.g., no need to redo symbolic factorization for
602: ! ILU/ICC preconditioners).
603: ! - If the nonzero structure of the matrix is different during
604: ! successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
605: ! must be used instead. If you are unsure whether the matrix
606: ! structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
607: ! - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
608: ! believes your assertion and does not check the structure
609: ! of the matrix. If you erroneously claim that the structure
610: ! is the same when it actually is not, the new preconditioner
611: ! will not function correctly. Thus, use this optimization
612: ! feature with caution!
614: flag = SAME_NONZERO_PATTERN
616: return
617: end
619: ! ---------------------------------------------------------------------
620: !
621: ! ApplicationJacobian - Computes Jacobian matrix, called by
622: ! the higher level routine FormJacobian().
623: !
624: ! Input Parameters:
625: ! x - local vector data
626: !
627: ! Output Parameters:
628: ! jac - Jacobian matrix
629: ! jac_prec - optionally different preconditioning matrix (not used here)
630: ! ierr - error code
631: !
632: ! Notes:
633: ! This routine uses standard Fortran-style computations over a 2-dim array.
634: !
635: subroutine ApplicationJacobian(x,jac,jac_prec,ierr)
636: implicit none
638: #include include/finclude/petsc.h
639: #include include/finclude/petscvec.h
640: #include include/finclude/petscmat.h
641: #include include/finclude/petscpc.h
642: #include include/finclude/petscsnes.h
644: ! Common blocks:
645: PetscReal lambda
646: PetscInt mx,my
647: common /params/ lambda,mx,my
649: ! Input/output variables:
650: PetscScalar x(mx,my)
651: Mat jac,jac_prec
652: PetscErrorCode ierr
654: ! Local variables:
655: PetscInt i,j,row(1),col(5),i1,i5
656: PetscScalar two,one, hx,hy
657: PetscScalar hxdhy,hydhx,sc,v(5)
659: ! Set parameters
661: i1 = 1
662: i5 = 5
663: one = 1.0
664: two = 2.0
665: hx = one/dble(mx-1)
666: hy = one/dble(my-1)
667: sc = hx*hy
668: hxdhy = hx/hy
669: hydhx = hy/hx
671: ! Compute entries of the Jacobian matrix
672: ! - Here, we set all entries for a particular row at once.
673: ! - Note that MatSetValues() uses 0-based row and column numbers
674: ! in Fortran as well as in C.
676: do 20 j=1,my
677: row(1) = (j-1)*mx - 1
678: do 10 i=1,mx
679: row(1) = row(1) + 1
680: ! boundary points
681: if (i .eq. 1 .or. j .eq. 1 &
682: & .or. i .eq. mx .or. j .eq. my) then
683: call MatSetValues(jac_prec,i1,row,i1,row,one, &
684: & INSERT_VALUES,ierr)
685: ! interior grid points
686: else
687: v(1) = -hxdhy
688: v(2) = -hydhx
689: v(3) = two*(hydhx + hxdhy) &
690: & - sc*lambda*exp(x(i,j))
691: v(4) = -hydhx
692: v(5) = -hxdhy
693: col(1) = row(1) - mx
694: col(2) = row(1) - 1
695: col(3) = row(1)
696: col(4) = row(1) + 1
697: col(5) = row(1) + mx
698: call MatSetValues(jac_prec,i1,row,i5,col,v, &
699: & INSERT_VALUES,ierr)
700: endif
701: 10 continue
702: 20 continue
704: return
705: end