Actual source code: ex5f90.F

petsc-3.4.5 2014-06-29
  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: !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 uniprocessor version of this code is snes/examples/tutorials/ex4f.F
 31: !
 32: !  --------------------------------------------------------------------------
 33: !  The following define must be used before including any PETSc include files
 34: !  into a module or interface. This is because they can't handle declarations
 35: !  in them
 36: !

 38:       module f90module
 39:       type userctx
 40: #include <finclude/petscsysdef.h>
 41: #include <finclude/petscvecdef.h>
 42: #include <finclude/petscdmdef.h>
 43:         PetscInt xs,xe,xm,gxs,gxe,gxm
 44:         PetscInt ys,ye,ym,gys,gye,gym
 45:         PetscInt mx,my
 46:         PetscMPIInt rank
 47:         double precision lambda
 48:       end type userctx

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

 75: #include <finclude/petscsys.h>
 76: #include <finclude/petscvec.h>
 77: #include <finclude/petscdmda.h>
 78: #include <finclude/petscis.h>
 79: #include <finclude/petscmat.h>
 80: #include <finclude/petscksp.h>
 81: #include <finclude/petscpc.h>
 82: #include <finclude/petscsnes.h>
 83: #include <finclude/petscvec.h90>
 84: #include <finclude/petscsnes.h90>

 86: !  Input/output variables:
 87:       SNES           snes
 88:       Vec            X,F
 89:       PetscErrorCode ierr
 90:       type (userctx) user
 91:       DM             da

 93: !  Declarations for use with local arrays:
 94:       PetscScalar,pointer :: lx_v(:),lf_v(:)
 95:       Vec            localX

 97: !  Scatter ghost points to local vector, using the 2-step process
 98: !     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
 99: !  By placing code between these two statements, computations can
100: !  be done while messages are in transition.
101:       call SNESGetDM(snes,da,ierr)
102:       call DMGetLocalVector(da,localX,ierr)
103:       call DMGlobalToLocalBegin(da,X,INSERT_VALUES,                     &
104:      &     localX,ierr)
105:       call DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)

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

115:       call VecGetArrayF90(localX,lx_v,ierr)
116:       call VecGetArrayF90(F,lf_v,ierr)

118: !  Compute function over the locally owned part of the grid
119:       call FormFunctionLocal(lx_v,lf_v,user,ierr)

121: !  Restore vectors
122:       call VecRestoreArrayF90(localX,lx_v,ierr)
123:       call VecRestoreArrayF90(F,lf_v,ierr)

125: !  Insert values into global vector

127:       call DMRestoreLocalVector(da,localX,ierr)
128:       call PetscLogFlops(11.0d0*user%ym*user%xm,ierr)

130: !      call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
131: !      call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)
132:       return
133:       end subroutine formfunction
134:       end module f90module

136:       module f90moduleinterfaces
137:         use f90module

139:       Interface SNESSetApplicationContext
140:         Subroutine SNESSetApplicationContext(snes,ctx,ierr)
141:         use f90module
142:           SNES snes
143:           type(userctx) ctx
144:           PetscErrorCode ierr
145:         End Subroutine
146:       End Interface SNESSetApplicationContext

148:       Interface SNESGetApplicationContext
149:         Subroutine SNESGetApplicationContext(snes,ctx,ierr)
150:         use f90module
151:           SNES snes
152:           type(userctx), pointer :: ctx
153:           PetscErrorCode ierr
154:         End Subroutine
155:       End Interface SNESGetApplicationContext
156:       end module f90moduleinterfaces

158:       program main
159:       use f90module
160:       use f90moduleinterfaces
161:       implicit none
162: !
163: #include <finclude/petscsys.h>
164: #include <finclude/petscvec.h>
165: #include <finclude/petscdmda.h>
166: #include <finclude/petscis.h>
167: #include <finclude/petscmat.h>
168: #include <finclude/petscksp.h>
169: #include <finclude/petscpc.h>
170: #include <finclude/petscsnes.h>
171: #include <finclude/petscvec.h90>
172: #include <finclude/petscdmda.h90>

174: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: !                   Variable declarations
176: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177: !
178: !  Variables:
179: !     snes        - nonlinear solver
180: !     x, r        - solution, residual vectors
181: !     J           - Jacobian matrix
182: !     its         - iterations for convergence
183: !     Nx, Ny      - number of preocessors in x- and y- directions
184: !     matrix_free - flag - 1 indicates matrix-free version
185: !
186:       SNES             snes
187:       Vec              x,r
188:       Mat              J
189:       PetscErrorCode   ierr
190:       PetscInt         its
191:       PetscBool        flg,matrix_free
192:       PetscInt         ione,nfour
193:       double precision lambda_max,lambda_min
194:       type (userctx)   user
195:       DM               da

197: !  Note: Any user-defined Fortran routines (such as FormJacobian)
198: !  MUST be declared as external.
199:       external FormInitialGuess,FormJacobian

201: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202: !  Initialize program
203: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
205:       call MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)

207: !  Initialize problem parameters
208:       lambda_max  = 6.81
209:       lambda_min  = 0.0
210:       user%lambda = 6.0
211:       ione = 1
212:       nfour = -4
213:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',             &
214:      &     user%lambda,flg,ierr)
215:       if (user%lambda .ge. lambda_max .or. user%lambda .le. lambda_min) &
216:      &     then
217:          if (user%rank .eq. 0) write(6,*) 'Lambda is out of range'
218:          SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
219:       endif

221: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222: !  Create nonlinear solver context
223: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

226: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227: !  Create vector data structures; set function evaluation routine
228: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

232: ! This really needs only the star-type stencil, but we use the box
233: ! stencil temporarily.
234:       call DMDACreate2d(PETSC_COMM_WORLD,                               &
235:      &     DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,                      &
236:      &     DMDA_STENCIL_BOX,nfour,nfour,PETSC_DECIDE,PETSC_DECIDE,          &
237:      &     ione,ione,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
238:       call DMDAGetInfo(da,PETSC_NULL_INTEGER,user%mx,user%my,           &
239:      &               PETSC_NULL_INTEGER,                                &
240:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
241:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
242:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
243:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
244:      &               PETSC_NULL_INTEGER,ierr)

246: !
247: !   Visualize the distribution of the array across the processors
248: !
249: !     call DMView(da,PETSC_VIEWER_DRAW_WORLD,ierr)

251: !  Extract global and local vectors from DMDA; then duplicate for remaining
252: !  vectors that are the same types
253:       call DMCreateGlobalVector(da,x,ierr)
254:       call VecDuplicate(x,r,ierr)

256: !  Get local grid boundaries (for 2-dimensional DMDA)
257:       call DMDAGetCorners(da,user%xs,user%ys,PETSC_NULL_INTEGER,        &
258:      &     user%xm,user%ym,PETSC_NULL_INTEGER,ierr)
259:       call DMDAGetGhostCorners(da,user%gxs,user%gys,                    &
260:      &     PETSC_NULL_INTEGER,user%gxm,user%gym,                        &
261:      &     PETSC_NULL_INTEGER,ierr)

263: !  Here we shift the starting indices up by one so that we can easily
264: !  use the Fortran convention of 1-based indices (rather 0-based indices).
265:       user%xs  = user%xs+1
266:       user%ys  = user%ys+1
267:       user%gxs = user%gxs+1
268:       user%gys = user%gys+1

270:       user%ye  = user%ys+user%ym-1
271:       user%xe  = user%xs+user%xm-1
272:       user%gye = user%gys+user%gym-1
273:       user%gxe = user%gxs+user%gxm-1

275:       call SNESSetApplicationContext(snes,user,ierr)

277: !  Set function evaluation routine and vector
278:       call SNESSetFunction(snes,r,FormFunction,user,ierr)

280: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
281: !  Create matrix data structure; set Jacobian evaluation routine
282: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

304:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_mf',         &
305:      &     matrix_free,ierr)
306:       if (.not. matrix_free) then
307:         call DMCreateMatrix(da,MATAIJ,J,ierr)
308:         call SNESSetJacobian(snes,J,J,FormJacobian,user,ierr)
309:       endif

311: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
312: !  Customize nonlinear solver; set runtime options
313: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
314: !  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
315:       call SNESSetDM(snes,da,ierr)
316:       call SNESSetFromOptions(snes,ierr)


319: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
320: !  Evaluate initial guess; then solve nonlinear system.
321: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
322: !  Note: The user should initialize the vector, x, with the initial guess
323: !  for the nonlinear solver prior to calling SNESSolve().  In particular,
324: !  to employ an initial guess of zero, the user should explicitly set
325: !  this vector to zero by calling VecSet().

327:       call FormInitialGuess(snes,x,ierr)
328:       call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
329:       call SNESGetIterationNumber(snes,its,ierr);
330:       if (user%rank .eq. 0) then
331:          write(6,100) its
332:       endif
333:   100 format('Number of SNES iterations = ',i5)

335: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
336: !  Free work space.  All PETSc objects should be destroyed when they
337: !  are no longer needed.
338: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
339:       if (.not. matrix_free) call MatDestroy(J,ierr)
340:       call VecDestroy(x,ierr)
341:       call VecDestroy(r,ierr)
342:       call SNESDestroy(snes,ierr)
343:       call DMDestroy(da,ierr)

345:       call PetscFinalize(ierr)
346:       end

348: ! ---------------------------------------------------------------------
349: !
350: !  FormInitialGuess - Forms initial approximation.
351: !
352: !  Input Parameters:
353: !  X - vector
354: !
355: !  Output Parameter:
356: !  X - vector
357: !
358: !  Notes:
359: !  This routine serves as a wrapper for the lower-level routine
360: !  "InitialGuessLocal", where the actual computations are
361: !  done using the standard Fortran style of treating the local
362: !  vector data as a multidimensional array over the local mesh.
363: !  This routine merely handles ghost point scatters and accesses
364: !  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
365: !
366:       subroutine FormInitialGuess(snes,X,ierr)
367:       use f90module
368:       use f90moduleinterfaces
369:       implicit none

371: #include <finclude/petscvec.h90>
372: #include <finclude/petscsys.h>
373: #include <finclude/petscvec.h>
374: #include <finclude/petscdmda.h>
375: #include <finclude/petscis.h>
376: #include <finclude/petscmat.h>
377: #include <finclude/petscksp.h>
378: #include <finclude/petscpc.h>
379: #include <finclude/petscsnes.h>

381: !  Input/output variables:
382:       SNES           snes
383:       type(userctx), pointer:: puser
384:       Vec            X
385:       PetscErrorCode ierr
386:       DM             da

388: !  Declarations for use with local arrays:
389:       PetscScalar,pointer :: lx_v(:)
390:       Vec               localX

392:       0
393:       call SNESGetDM(snes,da,ierr)
394:       call SNESGetApplicationContext(snes,puser,ierr)
395: !  Get a pointer to vector data.
396: !    - For default PETSc vectors, VecGetArray90() returns a pointer to
397: !      the data array. Otherwise, the routine is implementation dependent.
398: !    - You MUST call VecRestoreArrayF90() when you no longer need access to
399: !      the array.
400: !    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
401: !      and is useable from Fortran-90 Only.

403:       call DMGetLocalVector(da,localX,ierr)
404:       call VecGetArrayF90(localX,lx_v,ierr)

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

409: !  Restore vector
410:       call VecRestoreArrayF90(localX,lx_v,ierr)

412: !  Insert values into global vector
413:       call DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X,ierr)
414:       call DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X,ierr)
415:       call DMRestoreLocalVector(da,localX,ierr)

417:       return
418:       end

420: ! ---------------------------------------------------------------------
421: !
422: !  InitialGuessLocal - Computes initial approximation, called by
423: !  the higher level routine FormInitialGuess().
424: !
425: !  Input Parameter:
426: !  x - local vector data
427: !
428: !  Output Parameters:
429: !  x - local vector data
430: !  ierr - error code
431: !
432: !  Notes:
433: !  This routine uses standard Fortran-style computations over a 2-dim array.
434: !
435:       subroutine InitialGuessLocal(user,x,ierr)
436:       use f90module
437:       implicit none

439: #include <finclude/petscsys.h>
440: #include <finclude/petscvec.h>
441: #include <finclude/petscdmda.h>
442: #include <finclude/petscis.h>
443: #include <finclude/petscmat.h>
444: #include <finclude/petscksp.h>
445: #include <finclude/petscpc.h>
446: #include <finclude/petscsnes.h>

448: !  Input/output variables:
449:       type (userctx)         user
450:       PetscScalar  x(user%gxs:user%gxe,                                 &
451:      &              user%gys:user%gye)
452:       PetscErrorCode ierr

454: !  Local variables:
455:       PetscInt  i,j
456:       PetscScalar   temp1,temp,hx,hy
457:       PetscScalar   one

459: !  Set parameters

461:       0
462:       one    = 1.0
463:       hx     = one/(dble(user%mx-1))
464:       hy     = one/(dble(user%my-1))
465:       temp1  = user%lambda/(user%lambda + one)

467:       do 20 j=user%ys,user%ye
468:          temp = dble(min(j-1,user%my-j))*hy
469:          do 10 i=user%xs,user%xe
470:             if (i .eq. 1 .or. j .eq. 1                                  &
471:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
472:               x(i,j) = 0.0
473:             else
474:               x(i,j) = temp1 *                                          &
475:      &          sqrt(min(dble(min(i-1,user%mx-i)*hx),dble(temp)))
476:             endif
477:  10      continue
478:  20   continue

480:       return
481:       end

483: ! ---------------------------------------------------------------------
484: !
485: !  FormFunctionLocal - Computes nonlinear function, called by
486: !  the higher level routine FormFunction().
487: !
488: !  Input Parameter:
489: !  x - local vector data
490: !
491: !  Output Parameters:
492: !  f - local vector data, f(x)
493: !  ierr - error code
494: !
495: !  Notes:
496: !  This routine uses standard Fortran-style computations over a 2-dim array.
497: !
498:       subroutine FormFunctionLocal(x,f,user,ierr)
499:       use f90module

501:       implicit none

503: !  Input/output variables:
504:       type (userctx) user
505:       PetscScalar  x(user%gxs:user%gxe,                                         &
506:      &              user%gys:user%gye)
507:       PetscScalar  f(user%xs:user%xe,                                           &
508:      &              user%ys:user%ye)
509:       PetscErrorCode ierr

511: !  Local variables:
512:       PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
513:       PetscScalar u,uxx,uyy
514:       PetscInt  i,j

516:       one    = 1.0
517:       two    = 2.0
518:       hx     = one/dble(user%mx-1)
519:       hy     = one/dble(user%my-1)
520:       sc     = hx*hy*user%lambda
521:       hxdhy  = hx/hy
522:       hydhx  = hy/hx

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

526:       do 20 j=user%ys,user%ye
527:          do 10 i=user%xs,user%xe
528:             if (i .eq. 1 .or. j .eq. 1                                  &
529:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
530:                f(i,j) = x(i,j)
531:             else
532:                u = x(i,j)
533:                uxx = hydhx * (two*u                                     &
534:      &                - x(i-1,j) - x(i+1,j))
535:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
536:                f(i,j) = uxx + uyy - sc*exp(u)
537:             endif
538:  10      continue
539:  20   continue

541:       return
542:       end

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

595: #include <finclude/petscsys.h>
596: #include <finclude/petscvec.h>
597: #include <finclude/petscdmda.h>
598: #include <finclude/petscis.h>
599: #include <finclude/petscmat.h>
600: #include <finclude/petscksp.h>
601: #include <finclude/petscpc.h>
602: #include <finclude/petscsnes.h>

604: #include <finclude/petscvec.h90>

606: !  Input/output variables:
607:       SNES         snes
608:       Vec          X
609:       Mat          jac,jac_prec
610:       MatStructure flag
611:       type(userctx)  user
612:       PetscErrorCode ierr
613:       DM             da

615: !  Declarations for use with local arrays:
616:       PetscScalar,pointer :: lx_v(:)
617:       Vec            localX

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

624:       call SNESGetDM(snes,da,ierr)
625:       call DMGetLocalVector(da,localX,ierr)
626:       call DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,              &
627:      &     ierr)
628:       call DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)

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

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

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

641:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
642:       if (jac .ne. jac_prec) then
643:          call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
644:       endif
645:       call VecRestoreArrayF90(localX,lx_v,ierr)
646:       call DMRestoreLocalVector(da,localX,ierr)
647:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
648:       if (jac .ne. jac_prec) then
649:         call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
650:       endif

652: !  Set flag to indicate that the Jacobian matrix retains an identical
653: !  nonzero structure throughout all nonlinear iterations (although the
654: !  values of the entries change). Thus, we can save some work in setting
655: !  up the preconditioner (e.g., no need to redo symbolic factorization for
656: !  ILU/ICC preconditioners).
657: !   - If the nonzero structure of the matrix is different during
658: !     successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
659: !     must be used instead.  If you are unsure whether the matrix
660: !     structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
661: !   - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
662: !     believes your assertion and does not check the structure
663: !     of the matrix.  If you erroneously claim that the structure
664: !     is the same when it actually is not, the new preconditioner
665: !     will not function correctly.  Thus, use this optimization
666: !     feature with caution!

668:       flag = SAME_NONZERO_PATTERN

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

673:       call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,      &
674:      &                  ierr)

676:       return
677:       end

679: ! ---------------------------------------------------------------------
680: !
681: !  FormJacobianLocal - Computes Jacobian preconditioner matrix,
682: !  called by the higher level routine FormJacobian().
683: !
684: !  Input Parameters:
685: !  x        - local vector data
686: !
687: !  Output Parameters:
688: !  jac_prec - Jacobian preconditioner matrix
689: !  ierr     - error code
690: !
691: !  Notes:
692: !  This routine uses standard Fortran-style computations over a 2-dim array.
693: !
694: !  Notes:
695: !  Due to grid point reordering with DMDAs, we must always work
696: !  with the local grid points, and then transform them to the new
697: !  global numbering with the "ltog" mapping (via DMDAGetGlobalIndicesF90()).
698: !  We cannot work directly with the global numbers for the original
699: !  uniprocessor grid!
700: !
701: !  Two methods are available for imposing this transformation
702: !  when setting matrix entries:
703: !    (A) MatSetValuesLocal(), using the local ordering (including
704: !        ghost points!)
705: !        - Use DMDAGetGlobalIndicesF90() to extract the local-to-global map
706: !        - Associate this map with the matrix by calling
707: !          MatSetLocalToGlobalMapping() once
708: !        - Set matrix entries using the local ordering
709: !          by calling MatSetValuesLocal()
710: !    (B) MatSetValues(), using the global ordering
711: !        - Use DMDAGetGlobalIndicesF90() to extract the local-to-global map
712: !        - Then apply this map explicitly yourself
713: !        - Set matrix entries using the global ordering by calling
714: !          MatSetValues()
715: !  Option (A) seems cleaner/easier in many cases, and is the procedure
716: !  used in this example.
717: !
718:       subroutine FormJacobianLocal(x,jac_prec,user,ierr)
719:       use f90module
720:       implicit none

722: #include <finclude/petscsys.h>
723: #include <finclude/petscvec.h>
724: #include <finclude/petscdmda.h>
725: #include <finclude/petscis.h>
726: #include <finclude/petscmat.h>
727: #include <finclude/petscksp.h>
728: #include <finclude/petscpc.h>
729: #include <finclude/petscsnes.h>

731: !  Input/output variables:
732:       type (userctx) user
733:       PetscScalar    x(user%gxs:user%gxe,                                      &
734:      &               user%gys:user%gye)
735:       Mat            jac_prec
736:       PetscErrorCode ierr

738: !  Local variables:
739:       PetscInt    row,col(5),i,j
740:       PetscInt    ione,ifive
741:       PetscScalar two,one,hx,hy,hxdhy
742:       PetscScalar hydhx,sc,v(5)

744: !  Set parameters
745:       ione   = 1
746:       ifive  = 5
747:       one    = 1.0
748:       two    = 2.0
749:       hx     = one/dble(user%mx-1)
750:       hy     = one/dble(user%my-1)
751:       sc     = hx*hy
752:       hxdhy  = hx/hy
753:       hydhx  = hy/hx

755: !  Compute entries for the locally owned part of the Jacobian.
756: !   - Currently, all PETSc parallel matrix formats are partitioned by
757: !     contiguous chunks of rows across the processors.
758: !   - Each processor needs to insert only elements that it owns
759: !     locally (but any non-local elements will be sent to the
760: !     appropriate processor during matrix assembly).
761: !   - Here, we set all entries for a particular row at once.
762: !   - We can set matrix entries either using either
763: !     MatSetValuesLocal() or MatSetValues(), as discussed above.
764: !   - Note that MatSetValues() uses 0-based row and column numbers
765: !     in Fortran as well as in C.

767:       do 20 j=user%ys,user%ye
768:          row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
769:          do 10 i=user%xs,user%xe
770:             row = row + 1
771: !           boundary points
772:             if (i .eq. 1 .or. j .eq. 1                                  &
773:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
774:                col(1) = row
775:                v(1)   = one
776:                call MatSetValuesLocal(jac_prec,ione,row,ione,col,v,          &
777:      &                           INSERT_VALUES,ierr)
778: !           interior grid points
779:             else
780:                v(1) = -hxdhy
781:                v(2) = -hydhx
782:                v(3) = two*(hydhx + hxdhy)                               &
783:      &                  - sc*user%lambda*exp(x(i,j))
784:                v(4) = -hydhx
785:                v(5) = -hxdhy
786:                col(1) = row - user%gxm
787:                col(2) = row - 1
788:                col(3) = row
789:                col(4) = row + 1
790:                col(5) = row + user%gxm
791:                call MatSetValuesLocal(jac_prec,ione,row,ifive,col,v,         &
792:      &                                INSERT_VALUES,ierr)
793:             endif
794:  10      continue
795:  20   continue

797:       return
798:       end