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,ierr)

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

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

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

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

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

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

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

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

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

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

263: !  Set function evaluation routine and vector

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

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

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

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

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

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

304:       call SNESSetFromOptions(snes,ierr)

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

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

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

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

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

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

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

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

377:       0

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

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

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

392:       call InitialGuessLocal(user,lx_v,ierr)

394: !  Restore vector

396:       call VecRestoreArrayF90(localX,lx_v,ierr)

398: !  Insert values into global vector

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

403:       return
404:       end

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

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

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

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

444: !  Set parameters

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

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

468:       return
469:       end

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

489:       implicit none

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

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

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

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

512:       do 20 j=user%ys,user%ye
513:          do 10 i=user%xs,user%xe
514:             if (i .eq. 1 .or. j .eq. 1                                  &
515:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
516:                f(i,j) = x(i,j)
517:             else
518:                u = x(i,j)
519:                uxx = hydhx * (two*u                                     &
520:      &                - x(i-1,j) - x(i+1,j))
521:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
522:                f(i,j) = uxx + uyy - sc*exp(u)
523:             endif
524:  10      continue
525:  20   continue

527:       return
528:       end

530: ! ---------------------------------------------------------------------
531: !
532: !  FormJacobian - Evaluates Jacobian matrix.
533: !
534: !  Input Parameters:
535: !  snes     - the SNES context
536: !  x        - input vector
537: !  dummy    - optional user-defined context, as set by SNESSetJacobian()
538: !             (not used here)
539: !
540: !  Output Parameters:
541: !  jac      - Jacobian matrix
542: !  jac_prec - optionally different preconditioning matrix (not used here)
543: !  flag     - flag indicating matrix structure
544: !
545: !  Notes:
546: !  This routine serves as a wrapper for the lower-level routine
547: !  "FormJacobianLocal", where the actual computations are
548: !  done using the standard Fortran style of treating the local
549: !  vector data as a multidimensional array over the local mesh.
550: !  This routine merely accesses the local vector data via
551: !  VecGetArrayF90() and VecRestoreArrayF90().
552: !
553: !  Notes:
554: !  Due to grid point reordering with DAs, we must always work
555: !  with the local grid points, and then transform them to the new
556: !  global numbering with the "ltog" mapping (via DAGetGlobalIndicesF90()).
557: !  We cannot work directly with the global numbers for the original
558: !  uniprocessor grid!
559: !
560: !  Two methods are available for imposing this transformation
561: !  when setting matrix entries:
562: !    (A) MatSetValuesLocal(), using the local ordering (including
563: !        ghost points!)
564: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
565: !        - Associate this map with the matrix by calling
566: !          MatSetLocalToGlobalMapping() once
567: !        - Set matrix entries using the local ordering
568: !          by calling MatSetValuesLocal()
569: !    (B) MatSetValues(), using the global ordering
570: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
571: !        - Then apply this map explicitly yourself
572: !        - Set matrix entries using the global ordering by calling
573: !          MatSetValues()
574: !  Option (A) seems cleaner/easier in many cases, and is the procedure
575: !  used in this example.
576: !
577:       subroutine FormJacobian(snes,X,jac,jac_prec,flag,user,ierr)
578:       use f90module
579:       implicit none

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

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

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

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

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

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

615: !  Get a pointer to vector data

617:       call VecGetArrayF90(localX,lx_v,ierr)

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

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

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

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

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

649:       flag = SAME_NONZERO_PATTERN

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

654:        call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,ierr)

656:       return
657:       end

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

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

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

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

723: !  Set parameters

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

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

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

775:       return
776:       end