Actual source code: ex5f90.F

  1: ! "$Id: ex5f90.F,v 1.44 2001/09/11 18:47:20 bsmith Exp $";
  2: !
  3: !  Description: Solves a nonlinear system in parallel with SNES.
  4: !  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
  5: !  domain, using distributed arrays (DAs) to partition the parallel grid.
  6: !  The command line options include:
  7: !    -par <parameter>, where <parameter> indicates the nonlinearity of the problem
  8: !       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
  9: !
 10: !/*T
 11: !  Concepts: SNES^parallel Bratu example
 12: !  Concepts: DA^using distributed arrays;
 13: !  Processors: n
 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 uniprocessor version of this code is snes/examples/tutorials/ex4f.F
 32: !
 33: !  --------------------------------------------------------------------------
 34: !  The following define must be used before including any PETSc include files
 35: !  into a module or interface. This is because they can't handle declarations
 36: !  in them
 37: !

 39:       module f90module
 40:       type userctx
 41: #define PETSC_AVOID_DECLARATIONS
 42:  #include include/finclude/petsc.h
 43:  #include include/finclude/petscvec.h
 44:  #include include/finclude/petscda.h
 45: #undef PETSC_AVOID_DECLARATIONS
 46:         DA      da
 47:         integer xs,xe,xm,gxs,gxe,gxm
 48:         integer ys,ye,ym,gys,gye,gym
 49:         integer mx,my,rank
 50:         double precision lambda
 51:       end type userctx
 52:       contains
 53: ! ---------------------------------------------------------------------
 54: !
 55: !  FormFunction - Evaluates nonlinear function, F(x).
 56: !
 57: !  Input Parameters:
 58: !  snes - the SNES context
 59: !  X - input vector
 60: !  dummy - optional user-defined context, as set by SNESSetFunction()
 61: !          (not used here)
 62: !
 63: !  Output Parameter:
 64: !  F - function vector
 65: !
 66: !  Notes:
 67: !  This routine serves as a wrapper for the lower-level routine
 68: !  "FormFunctionLocal", where the actual computations are
 69: !  done using the standard Fortran style of treating the local
 70: !  vector data as a multidimensional array over the local mesh.
 71: !  This routine merely handles ghost point scatters and accesses
 72: !  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
 73: !
 74:       subroutine FormFunction(snes,X,F,user,ierr)
 75:       implicit none

 77:  #include include/finclude/petsc.h
 78:  #include include/finclude/petscvec.h
 79:  #include include/finclude/petscda.h
 80:  #include include/finclude/petscis.h
 81:  #include include/finclude/petscmat.h
 82:  #include include/finclude/petscksp.h
 83:  #include include/finclude/petscpc.h
 84:  #include include/finclude/petscsles.h
 85:  #include include/finclude/petscsnes.h

 87: #include "include/finclude/petscvec.h90"


 90: !  Input/output variables:
 91:       SNES           snes
 92:       Vec            X,F
 93:       integer        ierr
 94:       type (userctx) user

 96: !  Declarations for use with local arrays:
 97:       PetscScalar,pointer :: lx_v(:),lf_v(:)
 98:       Vec            localX

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

105:       call DAGetLocalVector(user%da,localX,ierr)
106:       call DAGlobalToLocalBegin(user%da,X,INSERT_VALUES,                &
107:      &     localX,ierr)
108:       call DAGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)

110: !  Get a pointer to vector data.
111: !    - For default PETSc vectors, VecGetArray90() returns a pointer to
112: !      the data array. Otherwise, the routine is implementation dependent.
113: !    - You MUST call VecRestoreArrayF90() when you no longer need access to
114: !      the array.
115: !    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
116: !      and is useable from Fortran-90 Only.

118:       call VecGetArrayF90(localX,lx_v,ierr)
119:       call VecGetArrayF90(F,lf_v,ierr)

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

123:       call FormFunctionLocal(lx_v,lf_v,user,ierr)

125: !  Restore vectors

127:       call VecRestoreArrayF90(localX,lx_v,ierr)
128:       call VecRestoreArrayF90(F,lf_v,ierr)

130: !  Insert values into global vector

132:       call DARestoreLocalVector(user%da,localX,ierr)
133:       call PetscLogFlops(11*user%ym*user%xm,ierr)

135: !      call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
136: !      call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)

138:       return
139:       end subroutine formfunction
140:       end module f90module



144:       program main
145:       use f90module
146:       implicit none
147: !
148: !
149:  #include include/finclude/petsc.h
150:  #include include/finclude/petscvec.h
151:  #include include/finclude/petscda.h
152:  #include include/finclude/petscis.h
153:  #include include/finclude/petscmat.h
154:  #include include/finclude/petscksp.h
155:  #include include/finclude/petscpc.h
156:  #include include/finclude/petscsles.h
157:  #include include/finclude/petscsnes.h
158: #include "include/finclude/petscvec.h90"
159: #include "include/finclude/petscda.h90"

161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162: !                   Variable declarations
163: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164: !
165: !  Variables:
166: !     snes        - nonlinear solver
167: !     x, r        - solution, residual vectors
168: !     J           - Jacobian matrix
169: !     its         - iterations for convergence
170: !     Nx, Ny      - number of preocessors in x- and y- directions
171: !     matrix_free - flag - 1 indicates matrix-free version
172: !
173: !
174:       SNES                   snes
175:       Vec                    x,r
176:       Mat                    J
177:       integer                its,matrix_free,flg,ierr
178:       double precision       lambda_max,lambda_min
179:       type (userctx)         user

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

184:       external FormInitialGuess,FormJacobian

186: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187: !  Initialize program
188: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

190:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
191:       call MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)

193: !  Initialize problem parameters

195:       lambda_max  = 6.81
196:       lambda_min  = 0.0
197:       user%lambda = 6.0
198:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',             &
199:      &     user%lambda,flg,ierr)
200:       if (user%lambda .ge. lambda_max .or. user%lambda .le. lambda_min) &
201:      &     then
202:          if (user%rank .eq. 0) write(6,*) 'Lambda is out of range'
203:          SETERRQ(1,' ',ierr)
204:       endif


207: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208: !  Create nonlinear solver context
209: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

211:       call SNESCreate(PETSC_COMM_WORLD,SNES_NONLINEAR_EQUATIONS,        &
212:      &                snes,ierr)

214: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215: !  Create vector data structures; set function evaluation routine
216: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

218: !  Create distributed array (DA) to manage parallel grid and vectors

220: ! This really needs only the star-type stencil, but we use the box
221: ! stencil temporarily.
222:       call DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,   &
223:      &     -4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,                         &
224:      &     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,user%da,ierr)
225:       call DAGetInfo(user%da,PETSC_NULL_INTEGER,user%mx,user%my,        &
226:      &               PETSC_NULL_INTEGER,                                &
227:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
228:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
229:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
230:      &               PETSC_NULL_INTEGER,ierr)
231: 
232: !
233: !   Visualize the distribution of the array across the processors
234: !
235: !     call DAView(user%da,PETSC_VIEWER_DRAW_WORLD,ierr)

237: !  Extract global and local vectors from DA; then duplicate for remaining
238: !  vectors that are the same types

240:       call DACreateGlobalVector(user%da,x,ierr)
241:       call VecDuplicate(x,r,ierr)

243: !  Get local grid boundaries (for 2-dimensional DA)

245:       call DAGetCorners(user%da,user%xs,user%ys,PETSC_NULL_INTEGER,     &
246:      &     user%xm,user%ym,PETSC_NULL_INTEGER,ierr)
247:       call DAGetGhostCorners(user%da,user%gxs,user%gys,                 &
248:      &     PETSC_NULL_INTEGER,user%gxm,user%gym,                        &
249:      &     PETSC_NULL_INTEGER,ierr)

251: !  Here we shift the starting indices up by one so that we can easily
252: !  use the Fortran convention of 1-based indices (rather 0-based indices).

254:       user%xs  = user%xs+1
255:       user%ys  = user%ys+1
256:       user%gxs = user%gxs+1
257:       user%gys = user%gys+1

259:       user%ye  = user%ys+user%ym-1
260:       user%xe  = user%xs+user%xm-1
261:       user%gye = user%gys+user%gym-1
262:       user%gxe = user%gxs+user%gxm-1

264: !  Set function evaluation routine and vector

266:       call SNESSetFunction(snes,r,FormFunction,user,ierr)

268: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269: !  Create matrix data structure; set Jacobian evaluation routine
270: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

272: !  Set Jacobian matrix data structure and default Jacobian evaluation
273: !  routine. User can override with:
274: !     -snes_fd : default finite differencing approximation of Jacobian
275: !     -snes_mf : matrix-free Newton-Krylov method with no preconditioning
276: !                (unless user explicitly sets preconditioner)
277: !     -snes_mf_operator : form preconditioning matrix as set by the user,
278: !                         but use matrix-free approx for Jacobian-vector
279: !                         products within Newton-Krylov method
280: !
281: !  Note:  For the parallel case, vectors and matrices MUST be partitioned
282: !     accordingly.  When using distributed arrays (DAs) to create vectors,
283: !     the DAs determine the problem partitioning.  We must explicitly
284: !     specify the local matrix dimensions upon its creation for compatibility
285: !     with the vector distribution.  Thus, the generic MatCreate() routine
286: !     is NOT sufficient when working with distributed arrays.
287: !
288: !     Note: Here we only approximately preallocate storage space for the
289: !     Jacobian.  See the users manual for a discussion of better techniques
290: !     for preallocating matrix memory.

292:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_mf',         &
293:      &     matrix_free,ierr)
294:       if (matrix_free .eq. 0) then
295:         call DAGetMatrix(user%da,MATMPIAIJ,J,ierr)
296:         call SNESSetJacobian(snes,J,J,FormJacobian,user,ierr)
297:       endif

299: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
300: !  Customize nonlinear solver; set runtime options
301: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

305:       call SNESSetFromOptions(snes,ierr)

307: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
308: !  Evaluate initial guess; then solve nonlinear system.
309: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

316:       call FormInitialGuess(user,x,ierr)
317:       call SNESSolve(snes,x,its,ierr)
318:       if (user%rank .eq. 0) then
319:          write(6,100) its
320:       endif
321:   100 format('Number of Newton iterations = ',i5)

323: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
324: !  Free work space.  All PETSc objects should be destroyed when they
325: !  are no longer needed.
326: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

328:       if (matrix_free .eq. 0) call MatDestroy(J,ierr)
329:       call VecDestroy(x,ierr)
330:       call VecDestroy(r,ierr)
331:       call SNESDestroy(snes,ierr)
332:       call DADestroy(user%da,ierr)
333:       call PetscFinalize(ierr)
334:       end

336: ! ---------------------------------------------------------------------
337: !
338: !  FormInitialGuess - Forms initial approximation.
339: !
340: !  Input Parameters:
341: !  X - vector
342: !
343: !  Output Parameter:
344: !  X - vector
345: !
346: !  Notes:
347: !  This routine serves as a wrapper for the lower-level routine
348: !  "InitialGuessLocal", where the actual computations are
349: !  done using the standard Fortran style of treating the local
350: !  vector data as a multidimensional array over the local mesh.
351: !  This routine merely handles ghost point scatters and accesses
352: !  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
353: !
354:       subroutine FormInitialGuess(user,X,ierr)
355:       use f90module
356:       implicit none

358: #include "include/finclude/petscvec.h90"
359:  #include include/finclude/petsc.h
360:  #include include/finclude/petscvec.h
361:  #include include/finclude/petscda.h
362:  #include include/finclude/petscis.h
363:  #include include/finclude/petscmat.h
364:  #include include/finclude/petscksp.h
365:  #include include/finclude/petscpc.h
366:  #include include/finclude/petscsles.h
367:  #include include/finclude/petscsnes.h

369: !  Input/output variables:
370:       type (userctx)         user
371:       Vec      X
372:       integer  ierr
373: 
374: !  Declarations for use with local arrays:
375:       PetscScalar,pointer :: lx_v(:)
376:       Vec               localX

378:       0

380: !  Get a pointer to vector data.
381: !    - For default PETSc vectors, VecGetArray90() returns a pointer to
382: !      the data array. Otherwise, the routine is implementation dependent.
383: !    - You MUST call VecRestoreArrayF90() when you no longer need access to
384: !      the array.
385: !    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
386: !      and is useable from Fortran-90 Only.

388:       call DAGetLocalVector(user%da,localX,ierr)
389:       call VecGetArrayF90(localX,lx_v,ierr)

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

393:       call InitialGuessLocal(user,lx_v,ierr)

395: !  Restore vector

397:       call VecRestoreArrayF90(localX,lx_v,ierr)

399: !  Insert values into global vector

401:       call DALocalToGlobal(user%da,localX,INSERT_VALUES,X,ierr)
402:       call DARestoreLocalVector(user%da,localX,ierr)

404:       return
405:       end

407: ! ---------------------------------------------------------------------
408: !
409: !  InitialGuessLocal - Computes initial approximation, called by
410: !  the higher level routine FormInitialGuess().
411: !
412: !  Input Parameter:
413: !  x - local vector data
414: !
415: !  Output Parameters:
416: !  x - local vector data
417: !  ierr - error code
418: !
419: !  Notes:
420: !  This routine uses standard Fortran-style computations over a 2-dim array.
421: !
422:       subroutine InitialGuessLocal(user,x,ierr)
423:       use f90module
424:       implicit none

426:  #include include/finclude/petsc.h
427:  #include include/finclude/petscvec.h
428:  #include include/finclude/petscda.h
429:  #include include/finclude/petscis.h
430:  #include include/finclude/petscmat.h
431:  #include include/finclude/petscksp.h
432:  #include include/finclude/petscpc.h
433:  #include include/finclude/petscsles.h
434:  #include include/finclude/petscsnes.h

436: !  Input/output variables:
437:       type (userctx)         user
438:       PetscScalar  x(user%gxs:user%gxe,user%gys:user%gye)
439:       integer ierr

441: !  Local variables:
442:       integer  i,j,hxdhy,hydhx
443:       PetscScalar   temp1,temp,hx,hy,sc,one

445: !  Set parameters

447:       ierr   = 0
448:       one    = 1.0
449:       hx     = one/(dble(user%mx-1))
450:       hy     = one/(dble(user%my-1))
451:       sc     = hx*hy*user%lambda
452:       hxdhy  = hx/hy
453:       hydhx  = hy/hx
454:       temp1  = user%lambda/(user%lambda + one)

456:       do 20 j=user%ys,user%ye
457:          temp = dble(min(j-1,user%my-j))*hy
458:          do 10 i=user%xs,user%xe
459:             if (i .eq. 1 .or. j .eq. 1                                  &
460:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
461:               x(i,j) = 0.0
462:             else
463:               x(i,j) = temp1 *                                          &
464:      &          sqrt(min(dble(min(i-1,user%mx-i)*hx),dble(temp)))
465:             endif
466:  10      continue
467:  20   continue

469:       return
470:       end

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

490:       implicit none

492: !  Input/output variables:
493:       type (userctx) user
494:       PetscScalar   x(user%gxs:user%gxe,user%gys:user%gye)
495:       PetscScalar   f(user%xs:user%xe,user%ys:user%ye)
496:       integer  ierr

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

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

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

513:       do 20 j=user%ys,user%ye
514:          do 10 i=user%xs,user%xe
515:             if (i .eq. 1 .or. j .eq. 1                                  &
516:      &             .or. i .eq. user%mx .or. j .eq. user%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: !  "FormJacobianLocal", 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: !  VecGetArrayF90() and VecRestoreArrayF90().
553: !
554: !  Notes:
555: !  Due to grid point reordering with DAs, we must always work
556: !  with the local grid points, and then transform them to the new
557: !  global numbering with the "ltog" mapping (via DAGetGlobalIndicesF90()).
558: !  We cannot work directly with the global numbers for the original
559: !  uniprocessor grid!
560: !
561: !  Two methods are available for imposing this transformation
562: !  when setting matrix entries:
563: !    (A) MatSetValuesLocal(), using the local ordering (including
564: !        ghost points!)
565: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
566: !        - Associate this map with the matrix by calling
567: !          MatSetLocalToGlobalMapping() once
568: !        - Set matrix entries using the local ordering
569: !          by calling MatSetValuesLocal()
570: !    (B) MatSetValues(), using the global ordering
571: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
572: !        - Then apply this map explicitly yourself
573: !        - Set matrix entries using the global ordering by calling
574: !          MatSetValues()
575: !  Option (A) seems cleaner/easier in many cases, and is the procedure
576: !  used in this example.
577: !
578:       subroutine FormJacobian(snes,X,jac,jac_prec,flag,user,ierr)
579:       use f90module
580:       implicit none

582:  #include include/finclude/petsc.h
583:  #include include/finclude/petscvec.h
584:  #include include/finclude/petscda.h
585:  #include include/finclude/petscis.h
586:  #include include/finclude/petscmat.h
587:  #include include/finclude/petscksp.h
588:  #include include/finclude/petscpc.h
589:  #include include/finclude/petscsles.h
590:  #include include/finclude/petscsnes.h

592: #include "include/finclude/petscvec.h90"

594: !  Input/output variables:
595:       SNES         snes
596:       Vec          X
597:       Mat          jac,jac_prec
598:       MatStructure flag
599:       type(userctx) user
600:       integer      ierr

602: !  Declarations for use with local arrays:
603:       PetscScalar,pointer :: lx_v(:)
604:       Vec            localX

606: !  Scatter ghost points to local vector, using the 2-step process
607: !     DAGlobalToLocalBegin(), DAGlobalToLocalEnd()
608: !  Computations can be done while messages are in transition,
609: !  by placing code between these two statements.

611:       call DAGetLocalVector(user%da,localX,ierr)
612:       call DAGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX,            &
613:      &     ierr)
614:       call DAGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)

616: !  Get a pointer to vector data

618:       call VecGetArrayF90(localX,lx_v,ierr)

620: !  Compute entries for the locally owned part of the Jacobian.

622:       call FormJacobianLocal(lx_v,jac,jac_prec,user,ierr)

624: !  Assemble matrix, using the 2-step process:
625: !     MatAssemblyBegin(), MatAssemblyEnd()
626: !  Computations can be done while messages are in transition,
627: !  by placing code between these two statements.

629:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
630:       call VecRestoreArrayF90(localX,lx_v,ierr)
631:       call DARestoreLocalVector(user%da,localX,ierr)
632:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)

634: !  Set flag to indicate that the Jacobian matrix retains an identical
635: !  nonzero structure throughout all nonlinear iterations (although the
636: !  values of the entries change). Thus, we can save some work in setting
637: !  up the preconditioner (e.g., no need to redo symbolic factorization for
638: !  ILU/ICC preconditioners).
639: !   - If the nonzero structure of the matrix is different during
640: !     successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
641: !     must be used instead.  If you are unsure whether the matrix
642: !     structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
643: !   - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
644: !     believes your assertion and does not check the structure
645: !     of the matrix.  If you erroneously claim that the structure
646: !     is the same when it actually is not, the new preconditioner
647: !     will not function correctly.  Thus, use this optimization
648: !     feature with caution!

650:       flag = SAME_NONZERO_PATTERN

652: !  Tell the matrix we will never add a new nonzero location to the
653: !  matrix. If we do it will generate an error.

655:        call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,ierr)

657:       return
658:       end

660: ! ---------------------------------------------------------------------
661: !
662: !  FormJacobianLocal - Computes Jacobian matrix, called by
663: !  the higher level routine FormJacobian().
664: !
665: !  Input Parameters:
666: !  x        - local vector data
667: !
668: !  Output Parameters:
669: !  jac      - Jacobian matrix
670: !  jac_prec - optionally different preconditioning matrix (not used here)
671: !  ierr     - error code
672: !
673: !  Notes:
674: !  This routine uses standard Fortran-style computations over a 2-dim array.
675: !
676: !  Notes:
677: !  Due to grid point reordering with DAs, we must always work
678: !  with the local grid points, and then transform them to the new
679: !  global numbering with the "ltog" mapping (via DAGetGlobalIndicesF90()).
680: !  We cannot work directly with the global numbers for the original
681: !  uniprocessor grid!
682: !
683: !  Two methods are available for imposing this transformation
684: !  when setting matrix entries:
685: !    (A) MatSetValuesLocal(), using the local ordering (including
686: !        ghost points!)
687: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
688: !        - Associate this map with the matrix by calling
689: !          MatSetLocalToGlobalMapping() once
690: !        - Set matrix entries using the local ordering
691: !          by calling MatSetValuesLocal()
692: !    (B) MatSetValues(), using the global ordering
693: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
694: !        - Then apply this map explicitly yourself
695: !        - Set matrix entries using the global ordering by calling
696: !          MatSetValues()
697: !  Option (A) seems cleaner/easier in many cases, and is the procedure
698: !  used in this example.
699: !
700:       subroutine FormJacobianLocal(x,jac,jac_prec,user,ierr)
701:       use f90module
702:       implicit none

704:  #include include/finclude/petsc.h
705:  #include include/finclude/petscvec.h
706:  #include include/finclude/petscda.h
707:  #include include/finclude/petscis.h
708:  #include include/finclude/petscmat.h
709:  #include include/finclude/petscksp.h
710:  #include include/finclude/petscpc.h
711:  #include include/finclude/petscsles.h
712:  #include include/finclude/petscsnes.h

714: !  Input/output variables:
715:       type (userctx) user
716:       PetscScalar   x(user%gxs:user%gxe,user%gys:user%gye)
717:       Mat      jac,jac_prec
718:       integer  ierr

720: !  Local variables:
721:       integer  row,col(5),i,j
722:       PetscScalar   two,one,hx,hy,hxdhy,hydhx,sc,v(5)

724: !  Set parameters

726:       one    = 1.0
727:       two    = 2.0
728:       hx     = one/dble(user%mx-1)
729:       hy     = one/dble(user%my-1)
730:       sc     = hx*hy
731:       hxdhy  = hx/hy
732:       hydhx  = hy/hx

734: !  Compute entries for the locally owned part of the Jacobian.
735: !   - Currently, all PETSc parallel matrix formats are partitioned by
736: !     contiguous chunks of rows across the processors.
737: !   - Each processor needs to insert only elements that it owns
738: !     locally (but any non-local elements will be sent to the
739: !     appropriate processor during matrix assembly).
740: !   - Here, we set all entries for a particular row at once.
741: !   - We can set matrix entries either using either
742: !     MatSetValuesLocal() or MatSetValues(), as discussed above.
743: !   - Note that MatSetValues() uses 0-based row and column numbers
744: !     in Fortran as well as in C.

746:       do 20 j=user%ys,user%ye
747:          row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
748:          do 10 i=user%xs,user%xe
749:             row = row + 1
750: !           boundary points
751:             if (i .eq. 1 .or. j .eq. 1                                  &
752:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
753:                call MatSetValuesLocal(jac,1,row,1,row,one,              &
754:      &                           INSERT_VALUES,ierr)
755: !           interior grid points
756:             else
757:                v(1) = -hxdhy
758:                v(2) = -hydhx
759:                v(3) = two*(hydhx + hxdhy)                               &
760:      &                  - sc*user%lambda*exp(x(i,j))
761:                v(4) = -hydhx
762:                v(5) = -hxdhy
763:                col(1) = row - user%gxm
764:                col(2) = row - 1
765:                col(3) = row
766:                col(4) = row + 1
767:                col(5) = row + user%gxm
768:                call MatSetValuesLocal(jac,1,row,5,col,v,                &
769:      &                                INSERT_VALUES,ierr)
770:             endif
771:  10      continue
772:  20   continue

774:       return
775:       end