Actual source code: ex5f90.F

  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:         DM      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:         double precision 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(snes,X,F,user,ierr)
 74:       implicit none

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

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

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

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

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

113:       call VecGetArrayF90(localX,lx_v,ierr)
114:       call VecGetArrayF90(F,lf_v,ierr)

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

119: !  Restore vectors
120:       call VecRestoreArrayF90(localX,lx_v,ierr)
121:       call VecRestoreArrayF90(F,lf_v,ierr)

123: !  Insert values into global vector

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

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

134:       module f90moduleinterfaces
135:         use f90module
136: 
137:       Interface SNESSetApplicationContext
138:         Subroutine SNESSetApplicationContext(snes,ctx,ierr)
139:         use f90module
140:           SNES snes
141:           type(userctx) ctx
142:           PetscErrorCode ierr
143:         End Subroutine
144:       End Interface SNESSetApplicationContext

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

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

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

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

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

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

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

224: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: !  Create vector data structures; set function evaluation routine
226: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

230: ! This really needs only the star-type stencil, but we use the box
231: ! stencil temporarily.
232:       call DMDACreate2d(PETSC_COMM_WORLD,                               &
233:      &     DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,                      &
234:      &     DMDA_STENCIL_BOX,nfour,nfour,PETSC_DECIDE,PETSC_DECIDE,          &
235:      &     ione,ione,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,user%da,ierr)
236:       call DMDAGetInfo(user%da,PETSC_NULL_INTEGER,user%mx,user%my,        &
237:      &               PETSC_NULL_INTEGER,                                &
238:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
239:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
240:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
241:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,             &
242:      &               PETSC_NULL_INTEGER,ierr)
243: 
244: !
245: !   Visualize the distribution of the array across the processors
246: !
247: !     call DMView(user%da,PETSC_VIEWER_DRAW_WORLD,ierr)

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

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

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

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

273:       call SNESSetApplicationContext(snes,user,ierr)

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

278: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
279: !  Create matrix data structure; set Jacobian evaluation routine
280: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

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

315: !     Test Fortran90 wrapper for SNESSet/Get ApplicationContext()
316:       call PetscOptionsGetBool(PETSC_NULL_CHARACTER,'-test_appctx',             &
317:      &     flg,PETSC_NULL_CHARACTER,ierr)
318:       if (flg) then
319:         call SNESGetApplicationContext(snes,puser,ierr)
320:         if (user%da .ne. puser%da) then
321:           write(*,*) "Error: uesr != puesr"
322:           write(*,*) "user: ", user
323:           write(*,*) "puesr: ", puser
324:         endif
325:       endif

327: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
328: !  Evaluate initial guess; then solve nonlinear system.
329: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
330: !  Note: The user should initialize the vector, x, with the initial guess
331: !  for the nonlinear solver prior to calling SNESSolve().  In particular,
332: !  to employ an initial guess of zero, the user should explicitly set
333: !  this vector to zero by calling VecSet().

335:       call FormInitialGuess(snes,x,ierr)
336:       call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
337:       call SNESGetIterationNumber(snes,its,ierr);
338:       if (user%rank .eq. 0) then
339:          write(6,100) its
340:       endif
341:   100 format('Number of Newton iterations = ',i5)

343: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
344: !  Free work space.  All PETSc objects should be destroyed when they
345: !  are no longer needed.
346: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
347:       if (.not. matrix_free) call MatDestroy(J,ierr)
348:       call VecDestroy(x,ierr)
349:       call VecDestroy(r,ierr)
350:       call SNESDestroy(snes,ierr)
351:       call DMDestroy(user%da,ierr)

353:       call PetscFinalize(ierr)
354:       end

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

379: #include <finclude/petscvec.h90>
380: #include <finclude/petscsys.h>
381: #include <finclude/petscvec.h>
382: #include <finclude/petscdmda.h>
383: #include <finclude/petscis.h>
384: #include <finclude/petscmat.h>
385: #include <finclude/petscksp.h>
386: #include <finclude/petscpc.h>
387: #include <finclude/petscsnes.h>

389: !  Input/output variables:
390:       SNES           snes
391:       type(userctx), pointer:: puser
392:       Vec            X
393:       PetscErrorCode ierr
394: 
395: !  Declarations for use with local arrays:
396:       PetscScalar,pointer :: lx_v(:)
397:       Vec               localX

399:       0
400:       call SNESGetApplicationContext(snes,puser,ierr)
401: !  Get a pointer to vector data.
402: !    - For default PETSc vectors, VecGetArray90() returns a pointer to
403: !      the data array. Otherwise, the routine is implementation dependent.
404: !    - You MUST call VecRestoreArrayF90() when you no longer need access to
405: !      the array.
406: !    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
407: !      and is useable from Fortran-90 Only.

409:       call DMGetLocalVector(puser%da,localX,ierr)
410:       call VecGetArrayF90(localX,lx_v,ierr)

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

415: !  Restore vector
416:       call VecRestoreArrayF90(localX,lx_v,ierr)

418: !  Insert values into global vector
419:       call DMLocalToGlobalBegin(puser%da,localX,INSERT_VALUES,X,ierr)
420:       call DMLocalToGlobalEnd(puser%da,localX,INSERT_VALUES,X,ierr)
421:       call DMRestoreLocalVector(puser%da,localX,ierr)

423:       return
424:       end

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

445: #include <finclude/petscsys.h>
446: #include <finclude/petscvec.h>
447: #include <finclude/petscdmda.h>
448: #include <finclude/petscis.h>
449: #include <finclude/petscmat.h>
450: #include <finclude/petscksp.h>
451: #include <finclude/petscpc.h>
452: #include <finclude/petscsnes.h>

454: !  Input/output variables:
455:       type (userctx)         user
456:       PetscScalar  x(user%gxs:user%gxe,                                 &
457:      &              user%gys:user%gye)
458:       PetscErrorCode ierr

460: !  Local variables:
461:       PetscInt  i,j
462:       PetscScalar   temp1,temp,hx,hy
463:       PetscScalar   one

465: !  Set parameters

467:       0
468:       one    = 1.0
469:       hx     = one/(dble(user%mx-1))
470:       hy     = one/(dble(user%my-1))
471:       temp1  = user%lambda/(user%lambda + one)

473:       do 20 j=user%ys,user%ye
474:          temp = dble(min(j-1,user%my-j))*hy
475:          do 10 i=user%xs,user%xe
476:             if (i .eq. 1 .or. j .eq. 1                                  &
477:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
478:               x(i,j) = 0.0
479:             else
480:               x(i,j) = temp1 *                                          &
481:      &          sqrt(min(dble(min(i-1,user%mx-i)*hx),dble(temp)))
482:             endif
483:  10      continue
484:  20   continue

486:       return
487:       end

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

507:       implicit none

509: !  Input/output variables:
510:       type (userctx) user
511:       PetscScalar  x(user%gxs:user%gxe,                                         &
512:      &              user%gys:user%gye)
513:       PetscScalar  f(user%xs:user%xe,                                           &
514:      &              user%ys:user%ye)
515:       PetscErrorCode ierr

517: !  Local variables:
518:       PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
519:       PetscScalar u,uxx,uyy
520:       PetscInt  i,j

522:       one    = 1.0
523:       two    = 2.0
524:       hx     = one/dble(user%mx-1)
525:       hy     = one/dble(user%my-1)
526:       sc     = hx*hy*user%lambda
527:       hxdhy  = hx/hy
528:       hydhx  = hy/hx

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

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

547:       return
548:       end

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

601: #include <finclude/petscsys.h>
602: #include <finclude/petscvec.h>
603: #include <finclude/petscdmda.h>
604: #include <finclude/petscis.h>
605: #include <finclude/petscmat.h>
606: #include <finclude/petscksp.h>
607: #include <finclude/petscpc.h>
608: #include <finclude/petscsnes.h>

610: #include <finclude/petscvec.h90>

612: !  Input/output variables:
613:       SNES         snes
614:       Vec          X
615:       Mat          jac,jac_prec
616:       MatStructure flag
617:       type(userctx)  user
618:       PetscErrorCode ierr

620: !  Declarations for use with local arrays:
621:       PetscScalar,pointer :: lx_v(:)
622:       Vec            localX

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

629:       call DMGetLocalVector(user%da,localX,ierr)
630:       call DMGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX,            &
631:      &     ierr)
632:       call DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)

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

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

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

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

672:       flag = SAME_NONZERO_PATTERN

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

677:       call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,      &
678:      &                  ierr)

680:       return
681:       end

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

726: #include <finclude/petscsys.h>
727: #include <finclude/petscvec.h>
728: #include <finclude/petscdmda.h>
729: #include <finclude/petscis.h>
730: #include <finclude/petscmat.h>
731: #include <finclude/petscksp.h>
732: #include <finclude/petscpc.h>
733: #include <finclude/petscsnes.h>

735: !  Input/output variables:
736:       type (userctx) user
737:       PetscScalar    x(user%gxs:user%gxe,                                      &
738:      &               user%gys:user%gye)
739:       Mat            jac_prec
740:       PetscErrorCode ierr

742: !  Local variables:
743:       PetscInt    row,col(5),i,j
744:       PetscInt    ione,ifive
745:       PetscScalar two,one,hx,hy,hxdhy
746:       PetscScalar hydhx,sc,v(5)

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

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

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

801:       return
802:       end