Actual source code: ex5f90t.F90

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: !
  2: !  Description: Solves a nonlinear system in parallel with SNES.
  3: !  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
  4: !  domain, using distributed arrays (DMDAs) to partition the parallel grid.
  5: !  The command line options include:
  6: !    -par <parameter>, where <parameter> indicates the nonlinearity of the problem
  7: !       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
  8: !
  9: !/*T
 10: !  Concepts: SNES^parallel Bratu example
 11: !  Concepts: DMDA^using distributed arrays;
 12: !  Processors: n
 13: !  TODO: Need to update to latest API
 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:  #include <petsc/finclude/petscdm.h>
 41:       use petscdmdef
 42:       type userctx
 43:         type(tDM) da
 44:         PetscInt xs,xe,xm,gxs,gxe,gxm
 45:         PetscInt ys,ye,ym,gys,gye,gym
 46:         PetscInt mx,my
 47:         PetscMPIInt rank
 48:         PetscReal lambda
 49:       end type userctx

 51:       contains
 52: ! ---------------------------------------------------------------------
 53: !
 54: !  FormFunction - Evaluates nonlinear function, F(x).
 55: !
 56: !  Input Parameters:
 57: !  snes - the SNES context
 58: !  X - input vector
 59: !  dummy - optional user-defined context, as set by SNESSetFunction()
 60: !          (not used here)
 61: !
 62: !  Output Parameter:
 63: !  F - function vector
 64: !
 65: !  Notes:
 66: !  This routine serves as a wrapper for the lower-level routine
 67: !  "FormFunctionLocal", where the actual computations are
 68: !  done using the standard Fortran style of treating the local
 69: !  vector data as a multidimensional array over the local mesh.
 70: !  This routine merely handles ghost point scatters and accesses
 71: !  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
 72: !
 73:       subroutine FormFunction(snesIn,X,F,user,ierr)
 74:  #include <petsc/finclude/petscsnes.h>
 75:       use petscsnes

 77: !  Input/output variables:
 78:       type(tSNES)     snesIn
 79:       type(tVec)      X,F
 80:       PetscErrorCode ierr
 81:       type (userctx) user

 83: !  Declarations for use with local arrays:
 84:       PetscScalar,pointer :: lx_v(:),lf_v(:)
 85:       type(tVec)              localX

 87: !  Scatter ghost points to local vector, using the 2-step process
 88: !     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
 89: !  By placing code between these two statements, computations can
 90: !  be done while messages are in transition.
 91:       call DMGetLocalVector(user%da,localX,ierr);CHKERRQ(ierr)
 92:       call DMGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX,ierr);CHKERRQ(ierr)
 93:       call DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr);CHKERRQ(ierr)

 95: !  Get a pointer to vector data.
 96: !    - For default PETSc vectors, VecGetArray90() returns a pointer to
 97: !      the data array. Otherwise, the routine is implementation dependent.
 98: !    - You MUST call VecRestoreArrayF90() when you no longer need access to
 99: !      the array.
100: !    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
101: !      and is useable from Fortran-90 Only.

103:       call VecGetArrayF90(localX,lx_v,ierr);CHKERRQ(ierr)
104:       call VecGetArrayF90(F,lf_v,ierr);CHKERRQ(ierr)

106: !  Compute function over the locally owned part of the grid
107:       call FormFunctionLocal(lx_v,lf_v,user,ierr);CHKERRQ(ierr)

109: !  Restore vectors
110:       call VecRestoreArrayF90(localX,lx_v,ierr);CHKERRQ(ierr)
111:       call VecRestoreArrayF90(F,lf_v,ierr);CHKERRQ(ierr)

113: !  Insert values into global vector

115:       call DMRestoreLocalVector(user%da,localX,ierr);CHKERRQ(ierr)
116:       call PetscLogFlops(11.0d0*user%ym*user%xm,ierr)

118: !      call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
119: !      call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)
120:       return
121:       end subroutine formfunction
122:       end module f90module

124:       module f90moduleinterfaces
125:         use f90module

127:       Interface SNESSetApplicationContext
128:         Subroutine SNESSetApplicationContext(snesIn,ctx,ierr)
129:  #include <petsc/finclude/petscsnes.h>
130:         use petscsnes
131:         use f90module
132:           type(tSNES)    snesIn
133:           type(userctx) ctx
134:           PetscErrorCode ierr
135:         End Subroutine
136:       End Interface SNESSetApplicationContext

138:       Interface SNESGetApplicationContext
139:         Subroutine SNESGetApplicationContext(snesIn,ctx,ierr)
140:  #include <petsc/finclude/petscsnes.h>
141:         use petscsnes
142:         use f90module
143:           type(tSNES)     snesIn
144:           type(userctx), pointer :: ctx
145:           PetscErrorCode ierr
146:         End Subroutine
147:       End Interface SNESGetApplicationContext
148:       end module f90moduleinterfaces

150:       program main
151:  #include <petsc/finclude/petscdm.h>
152:  #include <petsc/finclude/petscsnes.h>
153:       use petscdmda
154:       use petscdm
155:       use petscsnes
156:       use f90module
157:       use f90moduleinterfaces
158:       implicit none
159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: !                   Variable declarations
161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162: !
163: !  Variables:
164: !     mysnes      - nonlinear solver
165: !     x, r        - solution, residual vectors
166: !     J           - Jacobian matrix
167: !     its         - iterations for convergence
168: !     Nx, Ny      - number of preocessors in x- and y- directions
169: !     matrix_free - flag - 1 indicates matrix-free version
170: !
171:       type(tSNES)       mysnes
172:       type(tVec)        x,r
173:       type(tMat)        J
174:       PetscErrorCode   ierr
175:       PetscInt         its
176:       PetscBool        flg,matrix_free,set
177:       PetscInt         ione,nfour
178:       PetscReal lambda_max,lambda_min
179:       type(userctx)    user
180:       type(userctx), pointer:: puser
181:       type(tPetscOptions) :: options

183: !  Note: Any user-defined Fortran routines (such as FormJacobian)
184: !  MUST be declared as external.
185:       external FormInitialGuess,FormJacobian

187: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188: !  Initialize program
189: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
191:       if (ierr .ne. 0) then
192:         print*,'Unable to initialize PETSc'
193:         stop
194:       endif
195:       call MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)

197: !  Initialize problem parameters
198:       options%v = 0
199:       lambda_max  = 6.81
200:       lambda_min  = 0.0
201:       user%lambda = 6.0
202:       ione = 1
203:       nfour = 4
204:       call PetscOptionsGetReal(options,PETSC_NULL_CHARACTER,'-par',user%lambda,flg,ierr);CHKERRA(ierr)
205:       if (user%lambda .ge. lambda_max .or. user%lambda .le. lambda_min) then
206:          if (user%rank .eq. 0) write(6,*) 'Lambda is out of range'
207:          SETERRA(PETSC_COMM_SELF,1,' ')
208:       endif

210: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: !  Create nonlinear solver context
212: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213:       call SNESCreate(PETSC_COMM_WORLD,mysnes,ierr);CHKERRA(ierr)

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

219: !  Create distributed array (DMDA) to manage parallel grid and vectors

221: ! This really needs only the star-type stencil, but we use the box
222: ! stencil temporarily.
223:       call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nfour,nfour,PETSC_DECIDE,PETSC_DECIDE,ione,ione, &
224:      &     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,user%da,ierr);CHKERRA(ierr)
225:       call DMSetFromOptions(user%da,ierr);CHKERRA(ierr)
226:       call DMSetUp(user%da,ierr);CHKERRA(ierr)
227:       call DMDAGetInfo(user%da,PETSC_NULL_INTEGER,user%mx,user%my,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
228:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr)

230: !
231: !   Visualize the distribution of the array across the processors
232: !
233: !     call DMView(user%da,PETSC_VIEWER_DRAW_WORLD,ierr)

235: !  Extract global and local vectors from DMDA; then duplicate for remaining
236: !  vectors that are the same types
237:       call DMCreateGlobalVector(user%da,x,ierr);CHKERRA(ierr)
238:       call VecDuplicate(x,r,ierr);CHKERRA(ierr)

240: !  Get local grid boundaries (for 2-dimensional DMDA)
241:       call DMDAGetCorners(user%da,user%xs,user%ys,PETSC_NULL_INTEGER,user%xm,user%ym,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr)
242:       call DMDAGetGhostCorners(user%da,user%gxs,user%gys,PETSC_NULL_INTEGER,user%gxm,user%gym,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr)

244: !  Here we shift the starting indices up by one so that we can easily
245: !  use the Fortran convention of 1-based indices (rather 0-based indices).
246:       user%xs  = user%xs+1
247:       user%ys  = user%ys+1
248:       user%gxs = user%gxs+1
249:       user%gys = user%gys+1

251:       user%ye  = user%ys+user%ym-1
252:       user%xe  = user%xs+user%xm-1
253:       user%gye = user%gys+user%gym-1
254:       user%gxe = user%gxs+user%gxm-1

256:       call SNESSetApplicationContext(mysnes,user,ierr);CHKERRA(ierr)

258: !  Set function evaluation routine and vector
259:       call SNESSetFunction(mysnes,r,FormFunction,user,ierr);CHKERRA(ierr)

261: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262: !  Create matrix data structure; set Jacobian evaluation routine
263: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

285:       call PetscOptionsHasName(options,PETSC_NULL_CHARACTER,'-snes_mf',matrix_free,ierr);CHKERRA(ierr)
286:       if (.not. matrix_free) then
287:         call DMSetMatType(user%da,MATAIJ,ierr);CHKERRA(ierr)
288:         call DMCreateMatrix(user%da,J,ierr);CHKERRA(ierr)
289:         call SNESSetJacobian(mysnes,J,J,FormJacobian,user,ierr);CHKERRA(ierr)
290:       endif

292: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
293: !  Customize nonlinear solver; set runtime options
294: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
295: !  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
296:       call SNESSetFromOptions(mysnes,ierr);CHKERRA(ierr)

298: !     Test Fortran90 wrapper for SNESSet/Get ApplicationContext()
299:       call PetscOptionsGetBool(options,PETSC_NULL_CHARACTER,'-test_appctx',flg,set,ierr);CHKERRA(ierr)
300:       if (flg) then
301:         call SNESGetApplicationContext(mysnes,puser,ierr);CHKERRA(ierr)
302:       endif

304: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
305: !  Evaluate initial guess; then solve nonlinear system.
306: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307: !  Note: The user should initialize the vector, x, with the initial guess
308: !  for the nonlinear solver prior to calling SNESSolve().  In particular,
309: !  to employ an initial guess of zero, the user should explicitly set
310: !  this vector to zero by calling VecSet().

312:       call FormInitialGuess(mysnes,x,ierr);CHKERRA(ierr)
313:       call SNESSolve(mysnes,PETSC_NULL_VEC,x,ierr);CHKERRA(ierr)
314:       call SNESGetIterationNumber(mysnes,its,ierr);CHKERRA(ierr)
315:       if (user%rank .eq. 0) then
316:          write(6,100) its
317:       endif
318:   100 format('Number of SNES iterations = ',i5)

320: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
321: !  Free work space.  All PETSc objects should be destroyed when they
322: !  are no longer needed.
323: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
324:       if (.not. matrix_free) call MatDestroy(J,ierr);CHKERRA(ierr)
325:       call VecDestroy(x,ierr);CHKERRA(ierr)
326:       call VecDestroy(r,ierr);CHKERRA(ierr)
327:       call SNESDestroy(mysnes,ierr);CHKERRA(ierr)
328:       call DMDestroy(user%da,ierr);CHKERRA(ierr)

330:       call PetscFinalize(ierr)
331:       end

333: ! ---------------------------------------------------------------------
334: !
335: !  FormInitialGuess - Forms initial approximation.
336: !
337: !  Input Parameters:
338: !  X - vector
339: !
340: !  Output Parameter:
341: !  X - vector
342: !
343: !  Notes:
344: !  This routine serves as a wrapper for the lower-level routine
345: !  "InitialGuessLocal", where the actual computations are
346: !  done using the standard Fortran style of treating the local
347: !  vector data as a multidimensional array over the local mesh.
348: !  This routine merely handles ghost point scatters and accesses
349: !  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
350: !
351:       subroutine FormInitialGuess(mysnes,X,ierr)
352:  #include <petsc/finclude/petscsnes.h>
353:       use petscsnes
354:       use f90module
355:       use f90moduleinterfaces
356: !  Input/output variables:
357:       type(tSNES)     mysnes
358:       type(userctx), pointer:: puser
359:       type(tVec)      X
360:       PetscErrorCode ierr

362: !  Declarations for use with local arrays:
363:       PetscScalar,pointer :: lx_v(:)

365:       0
366:       call SNESGetApplicationContext(mysnes,puser,ierr)
367: !  Get a pointer to vector data.
368: !    - For default PETSc vectors, VecGetArray90() returns a pointer to
369: !      the data array. Otherwise, the routine is implementation dependent.
370: !    - You MUST call VecRestoreArrayF90() when you no longer need access to
371: !      the array.
372: !    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
373: !      and is useable from Fortran-90 Only.

375:       call VecGetArrayF90(X,lx_v,ierr)

377: !  Compute initial guess over the locally owned part of the grid
378:       call InitialGuessLocal(puser,lx_v,ierr)

380: !  Restore vector
381:       call VecRestoreArrayF90(X,lx_v,ierr)

383: !  Insert values into global vector

385:       return
386:       end

388: ! ---------------------------------------------------------------------
389: !
390: !  InitialGuessLocal - Computes initial approximation, called by
391: !  the higher level routine FormInitialGuess().
392: !
393: !  Input Parameter:
394: !  x - local vector data
395: !
396: !  Output Parameters:
397: !  x - local vector data
398: !  ierr - error code
399: !
400: !  Notes:
401: !  This routine uses standard Fortran-style computations over a 2-dim array.
402: !
403:       subroutine InitialGuessLocal(user,x,ierr)
404:  #include <petsc/finclude/petscsys.h>
405:       use petscsys
406:       use f90module
407: !  Input/output variables:
408:       type (userctx)         user
409:       PetscScalar  x(user%xs:user%xe,user%ys:user%ye)
410:       PetscErrorCode ierr

412: !  Local variables:
413:       PetscInt  i,j
414:       PetscScalar   temp1,temp,hx,hy
415:       PetscScalar   one

417: !  Set parameters

419:       0
420:       one    = 1.0
421:       hx     = one/(PetscIntToReal(user%mx-1))
422:       hy     = one/(PetscIntToReal(user%my-1))
423:       temp1  = user%lambda/(user%lambda + one)

425:       do 20 j=user%ys,user%ye
426:          temp = PetscIntToReal(min(j-1,user%my-j))*hy
427:          do 10 i=user%xs,user%xe
428:             if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
429:               x(i,j) = 0.0
430:             else
431:               x(i,j) = temp1 * sqrt(min(PetscIntToReal(min(i-1,user%mx-i)*hx),PetscIntToReal(temp)))
432:             endif
433:  10      continue
434:  20   continue

436:       return
437:       end

439: ! ---------------------------------------------------------------------
440: !
441: !  FormFunctionLocal - Computes nonlinear function, called by
442: !  the higher level routine FormFunction().
443: !
444: !  Input Parameter:
445: !  x - local vector data
446: !
447: !  Output Parameters:
448: !  f - local vector data, f(x)
449: !  ierr - error code
450: !
451: !  Notes:
452: !  This routine uses standard Fortran-style computations over a 2-dim array.
453: !
454:       subroutine FormFunctionLocal(x,f,user,ierr)
455:  #include <petsc/finclude/petscsys.h>
456:       use petscsys
457:       use f90module
458: !  Input/output variables:
459:       type (userctx) user
460:       PetscScalar  x(user%gxs:user%gxe,user%gys:user%gye)
461:       PetscScalar  f(user%xs:user%xe,user%ys:user%ye)
462:       PetscErrorCode ierr

464: !  Local variables:
465:       PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
466:       PetscScalar u,uxx,uyy
467:       PetscInt  i,j

469:       one    = 1.0
470:       two    = 2.0
471:       hx     = one/PetscIntToReal(user%mx-1)
472:       hy     = one/PetscIntToReal(user%my-1)
473:       sc     = hx*hy*user%lambda
474:       hxdhy  = hx/hy
475:       hydhx  = hy/hx

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

479:       do 20 j=user%ys,user%ye
480:          do 10 i=user%xs,user%xe
481:             if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
482:                f(i,j) = x(i,j)
483:             else
484:                u = x(i,j)
485:                uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j))
486:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
487:                f(i,j) = uxx + uyy - sc*exp(u)
488:             endif
489:  10      continue
490:  20   continue
491:       0
492:       return
493:       end

495: ! ---------------------------------------------------------------------
496: !
497: !  FormJacobian - Evaluates Jacobian matrix.
498: !
499: !  Input Parameters:
500: !  snes     - the SNES context
501: !  x        - input vector
502: !  dummy    - optional user-defined context, as set by SNESSetJacobian()
503: !             (not used here)
504: !
505: !  Output Parameters:
506: !  jac      - Jacobian matrix
507: !  jac_prec - optionally different preconditioning matrix (not used here)
508: !  flag     - flag indicating matrix structure
509: !
510: !  Notes:
511: !  This routine serves as a wrapper for the lower-level routine
512: !  "FormJacobianLocal", where the actual computations are
513: !  done using the standard Fortran style of treating the local
514: !  vector data as a multidimensional array over the local mesh.
515: !  This routine merely accesses the local vector data via
516: !  VecGetArrayF90() and VecRestoreArrayF90().
517: !
518: !  Notes:
519: !  Due to grid point reordering with DMDAs, we must always work
520: !  with the local grid points, and then transform them to the new
521: !  global numbering with the "ltog" mapping
522: !  We cannot work directly with the global numbers for the original
523: !  uniprocessor grid!
524: !
525: !  Two methods are available for imposing this transformation
526: !  when setting matrix entries:
527: !    (A) MatSetValuesLocal(), using the local ordering (including
528: !        ghost points!)
529: !        - Set matrix entries using the local ordering
530: !          by calling MatSetValuesLocal()
531: !    (B) MatSetValues(), using the global ordering
532: !        - Use DMGetLocalToGlobalMapping() then
533: !          ISLocalToGlobalMappingGetIndicesF90() to extract the local-to-global map
534: !        - Then apply this map explicitly yourself
535: !        - Set matrix entries using the global ordering by calling
536: !          MatSetValues()
537: !  Option (A) seems cleaner/easier in many cases, and is the procedure
538: !  used in this example.
539: !
540:       subroutine FormJacobian(mysnes,X,jac,jac_prec,user,ierr)
541:  #include <petsc/finclude/petscsnes.h>
542:       use petscsnes
543:       use f90module
544: !  Input/output variables:
545:       type(tSNES)     mysnes
546:       type(tVec)      X
547:       type(tMat)      jac,jac_prec
548:       type(userctx)  user
549:       PetscErrorCode ierr

551: !  Declarations for use with local arrays:
552:       PetscScalar,pointer :: lx_v(:)
553:       type(tVec)      localX

555: !  Scatter ghost points to local vector, using the 2-step process
556: !     DMGlobalToLocalBegin(), DMGlobalToLocalEnd()
557: !  Computations can be done while messages are in transition,
558: !  by placing code between these two statements.

560:       call DMGetLocalVector(user%da,localX,ierr)
561:       call DMGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX,ierr)
562:       call DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)

564: !  Get a pointer to vector data
565:       call VecGetArrayF90(localX,lx_v,ierr)

567: !  Compute entries for the locally owned part of the Jacobian preconditioner.
568:       call FormJacobianLocal(lx_v,jac_prec,user,ierr)

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

575:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
576: !      if (jac .ne. jac_prec) then
577:          call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
578: !      endif
579:       call VecRestoreArrayF90(localX,lx_v,ierr)
580:       call DMRestoreLocalVector(user%da,localX,ierr)
581:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
582: !      if (jac .ne. jac_prec) then
583:         call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
584: !      endif

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

589:       call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,ierr)

591:       return
592:       end

594: ! ---------------------------------------------------------------------
595: !
596: !  FormJacobianLocal - Computes Jacobian preconditioner matrix,
597: !  called by the higher level routine FormJacobian().
598: !
599: !  Input Parameters:
600: !  x        - local vector data
601: !
602: !  Output Parameters:
603: !  jac_prec - Jacobian preconditioner matrix
604: !  ierr     - error code
605: !
606: !  Notes:
607: !  This routine uses standard Fortran-style computations over a 2-dim array.
608: !
609: !  Notes:
610: !  Due to grid point reordering with DMDAs, we must always work
611: !  with the local grid points, and then transform them to the new
612: !  global numbering with the "ltog" mapping
613: !  We cannot work directly with the global numbers for the original
614: !  uniprocessor grid!
615: !
616: !  Two methods are available for imposing this transformation
617: !  when setting matrix entries:
618: !    (A) MatSetValuesLocal(), using the local ordering (including
619: !        ghost points!)
620: !        - Set matrix entries using the local ordering
621: !          by calling MatSetValuesLocal()
622: !    (B) MatSetValues(), using the global ordering
623: !        - Set matrix entries using the global ordering by calling
624: !          MatSetValues()
625: !  Option (A) seems cleaner/easier in many cases, and is the procedure
626: !  used in this example.
627: !
628:       subroutine FormJacobianLocal(x,jac_prec,user,ierr)
629:  #include <petsc/finclude/petscmat.h>
630:       use petscmat
631:       use f90module
632: !  Input/output variables:
633:       type (userctx) user
634:       PetscScalar    x(user%gxs:user%gxe,user%gys:user%gye)
635:       type(tMat)      jac_prec
636:       PetscErrorCode ierr

638: !  Local variables:
639:       PetscInt    row,col(5),i,j
640:       PetscInt    ione,ifive
641:       PetscScalar two,one,hx,hy,hxdhy
642:       PetscScalar hydhx,sc,v(5)

644: !  Set parameters
645:       ione   = 1
646:       ifive  = 5
647:       one    = 1.0
648:       two    = 2.0
649:       hx     = one/PetscIntToReal(user%mx-1)
650:       hy     = one/PetscIntToReal(user%my-1)
651:       sc     = hx*hy
652:       hxdhy  = hx/hy
653:       hydhx  = hy/hx

655: !  Compute entries for the locally owned part of the Jacobian.
656: !   - Currently, all PETSc parallel matrix formats are partitioned by
657: !     contiguous chunks of rows across the processors.
658: !   - Each processor needs to insert only elements that it owns
659: !     locally (but any non-local elements will be sent to the
660: !     appropriate processor during matrix assembly).
661: !   - Here, we set all entries for a particular row at once.
662: !   - We can set matrix entries either using either
663: !     MatSetValuesLocal() or MatSetValues(), as discussed above.
664: !   - Note that MatSetValues() uses 0-based row and column numbers
665: !     in Fortran as well as in C.

667:       do 20 j=user%ys,user%ye
668:          row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
669:          do 10 i=user%xs,user%xe
670:             row = row + 1
671: !           boundary points
672:             if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
673:                col(1) = row
674:                v(1)   = one
675:                call MatSetValuesLocal(jac_prec,ione,row,ione,col,v,INSERT_VALUES,ierr)
676: !           interior grid points
677:             else
678:                v(1) = -hxdhy
679:                v(2) = -hydhx
680:                v(3) = two*(hydhx + hxdhy) - sc*user%lambda*exp(x(i,j))
681:                v(4) = -hydhx
682:                v(5) = -hxdhy
683:                col(1) = row - user%gxm
684:                col(2) = row - 1
685:                col(3) = row
686:                col(4) = row + 1
687:                col(5) = row + user%gxm
688:                call MatSetValuesLocal(jac_prec,ione,row,ifive,col,v,INSERT_VALUES,ierr)
689:             endif
690:  10      continue
691:  20   continue
692:       return
693:       end