Actual source code: ex1f.F
petsc-3.4.2 2013-07-02
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 <finclude/petscsys.h>
57: #include <finclude/petscvec.h>
58: #include <finclude/petscis.h>
59: #include <finclude/petscdraw.h>
60: #include <finclude/petscmat.h>
61: #include <finclude/petscksp.h>
62: #include <finclude/petscpc.h>
63: #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: PetscBool 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
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_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(PETSC_COMM_SELF,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 (.not. matrix_free) 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 = .false.
179: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_fd_coloring', &
180: & fd_coloring,ierr)
181: if (fd_coloring) then
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,MATCOLORINGNATURAL,iscoloring,ierr)
195: !
196: ! Create the data structure that SNESComputeJacobianDefaultColor() uses
197: ! to compute the actual Jacobians via finite differences.
198: !
199: call MatFDColoringCreate(J,iscoloring,fdcoloring,ierr)
200: call MatFDColoringSetFunction(fdcoloring,FormFunction, &
201: & PETSC_NULL_OBJECT,ierr)
202: call MatFDColoringSetFromOptions(fdcoloring,ierr)
203: !
204: ! Tell SNES to use the routine SNESComputeJacobianDefaultColor()
205: ! to compute Jacobians.
206: !
207: call SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor, &
208: & fdcoloring,ierr)
209: call ISColoringDestroy(iscoloring,ierr)
211: else if (.not. matrix_free) 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 SNES iterations = ',i1)
256: ! PetscDraw contour plot of solution
258: call PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER, &
259: & 'Solution',300,0,300,300,draw,ierr)
260: call PetscDrawSetFromOptions(draw,ierr)
262: call VecGetArray(x,lx_v,lx_i,ierr)
263: call PetscDrawTensorContour(draw,mx,my,PETSC_NULL_DOUBLE, &
264: & PETSC_NULL_DOUBLE,lx_v(lx_i+1),ierr)
265: call VecRestoreArray(x,lx_v,lx_i,ierr)
267: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268: ! Free work space. All PETSc objects should be destroyed when they
269: ! are no longer needed.
270: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272: if (.not. matrix_free) call MatDestroy(J,ierr)
273: if (fd_coloring) call MatFDColoringDestroy(fdcoloring,ierr)
275: call VecDestroy(x,ierr)
276: call VecDestroy(r,ierr)
277: call SNESDestroy(snes,ierr)
278: call PetscDrawDestroy(draw,ierr)
279: call PetscFinalize(ierr)
280: end
282: ! ---------------------------------------------------------------------
283: !
284: ! FormInitialGuess - Forms initial approximation.
285: !
286: ! Input Parameter:
287: ! X - vector
288: !
289: ! Output Parameters:
290: ! X - vector
291: ! ierr - error code
292: !
293: ! Notes:
294: ! This routine serves as a wrapper for the lower-level routine
295: ! "ApplicationInitialGuess", where the actual computations are
296: ! done using the standard Fortran style of treating the local
297: ! vector data as a multidimensional array over the local mesh.
298: ! This routine merely accesses the local vector data via
299: ! VecGetArray() and VecRestoreArray().
300: !
301: subroutine FormInitialGuess(X,ierr)
302: implicit none
304: #include <finclude/petscsys.h>
305: #include <finclude/petscvec.h>
306: #include <finclude/petscmat.h>
307: #include <finclude/petscsnes.h>
309: ! Input/output variables:
310: Vec X
311: PetscErrorCode ierr
313: ! Declarations for use with local arrays:
314: PetscScalar lx_v(0:1)
315: PetscOffset lx_i
317: 0
319: ! Get a pointer to vector data.
320: ! - For default PETSc vectors, VecGetArray() returns a pointer to
321: ! the data array. Otherwise, the routine is implementation dependent.
322: ! - You MUST call VecRestoreArray() when you no longer need access to
323: ! the array.
324: ! - Note that the Fortran interface to VecGetArray() differs from the
325: ! C version. See the users manual for details.
327: call VecGetArray(X,lx_v,lx_i,ierr)
329: ! Compute initial guess
331: call ApplicationInitialGuess(lx_v(lx_i),ierr)
333: ! Restore vector
335: call VecRestoreArray(X,lx_v,lx_i,ierr)
337: return
338: end
340: ! ---------------------------------------------------------------------
341: !
342: ! ApplicationInitialGuess - Computes initial approximation, called by
343: ! the higher level routine FormInitialGuess().
344: !
345: ! Input Parameter:
346: ! x - local vector data
347: !
348: ! Output Parameters:
349: ! f - local vector data, f(x)
350: ! ierr - error code
351: !
352: ! Notes:
353: ! This routine uses standard Fortran-style computations over a 2-dim array.
354: !
355: subroutine ApplicationInitialGuess(x,ierr)
357: implicit none
359: ! Common blocks:
360: PetscReal lambda
361: PetscInt mx,my
362: common /params/ lambda,mx,my
364: ! Input/output variables:
365: PetscScalar x(mx,my)
366: PetscErrorCode ierr
368: ! Local variables:
369: PetscInt i,j
370: PetscScalar temp1,temp,hx,hy,one
372: ! Set parameters
374: 0
375: one = 1.0
376: hx = one/(dble(mx-1))
377: hy = one/(dble(my-1))
378: temp1 = lambda/(lambda + one)
380: do 20 j=1,my
381: temp = dble(min(j-1,my-j))*hy
382: do 10 i=1,mx
383: if (i .eq. 1 .or. j .eq. 1 &
384: & .or. i .eq. mx .or. j .eq. my) then
385: x(i,j) = 0.0
386: else
387: x(i,j) = temp1 * &
388: & sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
389: endif
390: 10 continue
391: 20 continue
393: return
394: end
396: ! ---------------------------------------------------------------------
397: !
398: ! FormFunction - Evaluates nonlinear function, F(x).
399: !
400: ! Input Parameters:
401: ! snes - the SNES context
402: ! X - input vector
403: ! dummy - optional user-defined context, as set by SNESSetFunction()
404: ! (not used here)
405: !
406: ! Output Parameter:
407: ! F - vector with newly computed function
408: !
409: ! Notes:
410: ! This routine serves as a wrapper for the lower-level routine
411: ! "ApplicationFunction", where the actual computations are
412: ! done using the standard Fortran style of treating the local
413: ! vector data as a multidimensional array over the local mesh.
414: ! This routine merely accesses the local vector data via
415: ! VecGetArray() and VecRestoreArray().
416: !
417: subroutine FormFunction(snes,X,F,dummy,ierr)
418: implicit none
420: #include <finclude/petscsys.h>
421: #include <finclude/petscvec.h>
422: #include <finclude/petscsnes.h>
424: ! Input/output variables:
425: SNES snes
426: Vec X,F
427: PetscFortranAddr dummy
428: PetscErrorCode ierr
430: ! Common blocks:
431: PetscReal lambda
432: PetscInt mx,my
433: common /params/ lambda,mx,my
435: ! Declarations for use with local arrays:
436: PetscScalar lx_v(0:1),lf_v(0:1)
437: PetscOffset lx_i,lf_i
439: ! Get pointers to vector data.
440: ! - For default PETSc vectors, VecGetArray() returns a pointer to
441: ! the data array. Otherwise, the routine is implementation dependent.
442: ! - You MUST call VecRestoreArray() when you no longer need access to
443: ! the array.
444: ! - Note that the Fortran interface to VecGetArray() differs from the
445: ! C version. See the Fortran chapter of the users manual for details.
447: call VecGetArray(X,lx_v,lx_i,ierr)
448: call VecGetArray(F,lf_v,lf_i,ierr)
450: ! Compute function
452: call ApplicationFunction(lx_v(lx_i),lf_v(lf_i),ierr)
454: ! Restore vectors
456: call VecRestoreArray(X,lx_v,lx_i,ierr)
457: call VecRestoreArray(F,lf_v,lf_i,ierr)
459: call PetscLogFlops(11.0d0*mx*my,ierr)
461: return
462: end
464: ! ---------------------------------------------------------------------
465: !
466: ! ApplicationFunction - Computes nonlinear function, called by
467: ! the higher level routine FormFunction().
468: !
469: ! Input Parameter:
470: ! x - local vector data
471: !
472: ! Output Parameters:
473: ! f - local vector data, f(x)
474: ! ierr - error code
475: !
476: ! Notes:
477: ! This routine uses standard Fortran-style computations over a 2-dim array.
478: !
479: subroutine ApplicationFunction(x,f,ierr)
481: implicit none
483: ! Common blocks:
484: PetscReal lambda
485: PetscInt mx,my
486: common /params/ lambda,mx,my
488: ! Input/output variables:
489: PetscScalar x(mx,my),f(mx,my)
490: PetscErrorCode ierr
492: ! Local variables:
493: PetscScalar two,one,hx,hy
494: PetscScalar hxdhy,hydhx,sc
495: PetscScalar u,uxx,uyy
496: PetscInt i,j
498: 0
499: one = 1.0
500: two = 2.0
501: hx = one/dble(mx-1)
502: hy = one/dble(my-1)
503: sc = hx*hy*lambda
504: hxdhy = hx/hy
505: hydhx = hy/hx
507: ! Compute function
509: do 20 j=1,my
510: do 10 i=1,mx
511: if (i .eq. 1 .or. j .eq. 1 &
512: & .or. i .eq. mx .or. j .eq. my) then
513: f(i,j) = x(i,j)
514: else
515: u = x(i,j)
516: uxx = hydhx * (two*u &
517: & - x(i-1,j) - x(i+1,j))
518: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
519: f(i,j) = uxx + uyy - sc*exp(u)
520: endif
521: 10 continue
522: 20 continue
524: return
525: end
527: ! ---------------------------------------------------------------------
528: !
529: ! FormJacobian - Evaluates Jacobian matrix.
530: !
531: ! Input Parameters:
532: ! snes - the SNES context
533: ! x - input vector
534: ! dummy - optional user-defined context, as set by SNESSetJacobian()
535: ! (not used here)
536: !
537: ! Output Parameters:
538: ! jac - Jacobian matrix
539: ! jac_prec - optionally different preconditioning matrix (not used here)
540: ! flag - flag indicating matrix structure
541: !
542: ! Notes:
543: ! This routine serves as a wrapper for the lower-level routine
544: ! "ApplicationJacobian", where the actual computations are
545: ! done using the standard Fortran style of treating the local
546: ! vector data as a multidimensional array over the local mesh.
547: ! This routine merely accesses the local vector data via
548: ! VecGetArray() and VecRestoreArray().
549: !
550: subroutine FormJacobian(snes,X,jac,jac_prec,flag,dummy,ierr)
551: implicit none
553: #include <finclude/petscsys.h>
554: #include <finclude/petscvec.h>
555: #include <finclude/petscmat.h>
556: #include <finclude/petscpc.h>
557: #include <finclude/petscsnes.h>
559: ! Input/output variables:
560: SNES snes
561: Vec X
562: Mat jac,jac_prec
563: MatStructure flag
564: PetscErrorCode ierr
565: integer dummy
567: ! Common blocks:
568: PetscReal lambda
569: PetscInt mx,my
570: common /params/ lambda,mx,my
572: ! Declarations for use with local array:
573: PetscScalar lx_v(0:1)
574: PetscOffset lx_i
576: ! Get a pointer to vector data
578: call VecGetArray(X,lx_v,lx_i,ierr)
580: ! Compute Jacobian entries
582: call ApplicationJacobian(lx_v(lx_i),jac,jac_prec,ierr)
584: ! Restore vector
586: call VecRestoreArray(X,lx_v,lx_i,ierr)
588: ! Assemble matrix
590: call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
591: call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
593: ! Set flag to indicate that the Jacobian matrix retains an identical
594: ! nonzero structure throughout all nonlinear iterations (although the
595: ! values of the entries change). Thus, we can save some work in setting
596: ! up the preconditioner (e.g., no need to redo symbolic factorization for
597: ! ILU/ICC preconditioners).
598: ! - If the nonzero structure of the matrix is different during
599: ! successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
600: ! must be used instead. If you are unsure whether the matrix
601: ! structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
602: ! - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
603: ! believes your assertion and does not check the structure
604: ! of the matrix. If you erroneously claim that the structure
605: ! is the same when it actually is not, the new preconditioner
606: ! will not function correctly. Thus, use this optimization
607: ! feature with caution!
609: flag = SAME_NONZERO_PATTERN
611: return
612: end
614: ! ---------------------------------------------------------------------
615: !
616: ! ApplicationJacobian - Computes Jacobian matrix, called by
617: ! the higher level routine FormJacobian().
618: !
619: ! Input Parameters:
620: ! x - local vector data
621: !
622: ! Output Parameters:
623: ! jac - Jacobian matrix
624: ! jac_prec - optionally different preconditioning matrix (not used here)
625: ! ierr - error code
626: !
627: ! Notes:
628: ! This routine uses standard Fortran-style computations over a 2-dim array.
629: !
630: subroutine ApplicationJacobian(x,jac,jac_prec,ierr)
631: implicit none
633: #include <finclude/petscsys.h>
634: #include <finclude/petscvec.h>
635: #include <finclude/petscmat.h>
636: #include <finclude/petscpc.h>
637: #include <finclude/petscsnes.h>
639: ! Common blocks:
640: PetscReal lambda
641: PetscInt mx,my
642: common /params/ lambda,mx,my
644: ! Input/output variables:
645: PetscScalar x(mx,my)
646: Mat jac,jac_prec
647: PetscErrorCode ierr
649: ! Local variables:
650: PetscInt i,j,row(1),col(5),i1,i5
651: PetscScalar two,one, hx,hy
652: PetscScalar hxdhy,hydhx,sc,v(5)
654: ! Set parameters
656: i1 = 1
657: i5 = 5
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(1) = (j-1)*mx - 1
673: do 10 i=1,mx
674: row(1) = row(1) + 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,i1,row,i1,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(1) - mx
689: col(2) = row(1) - 1
690: col(3) = row(1)
691: col(4) = row(1) + 1
692: col(5) = row(1) + mx
693: call MatSetValues(jac_prec,i1,row,i5,col,v, &
694: & INSERT_VALUES,ierr)
695: endif
696: 10 continue
697: 20 continue
699: return
700: end