Actual source code: ex39f90.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 (DAs) 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: !   Modified from ex5f90.F by Mike McCourt <mccomic@iit.edu>
 10: !   for testing Fortran interface on
 11: !   SNESLineSearchSet(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
 12: !
 13: !  --------------------------------------------------------------------------
 14: !
 15: !  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 16: !  the partial differential equation
 17: !
 18: !          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 19: !
 20: !  with boundary conditions
 21: !
 22: !           u = 0  for  x = 0, x = 1, y = 0, y = 1.
 23: !
 24: !  A finite difference approximation with the usual 5-point stencil
 25: !  is used to discretize the boundary value problem to obtain a nonlinear
 26: !  system of equations.
 27: !
 28: !  The uniprocessor version of this code is snes/examples/tutorials/ex4f.F
 29: !
 30: !  --------------------------------------------------------------------------
 31: !  The following define must be used before including any PETSc include files
 32: !  into a module or interface. This is because they can't handle declarations
 33: !  in them
 34: !

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

 72:  #include finclude/petscsys.h
 73:  #include finclude/petscvec.h
 74:  #include finclude/petscda.h
 75:  #include finclude/petscis.h
 76:  #include finclude/petscmat.h
 77:  #include finclude/petscksp.h
 78:  #include finclude/petscpc.h
 79:  #include finclude/petscsnes.h

 81: #include "finclude/petscvec.h90"


 84: !  Input/output variables:
 85:       SNES           snes
 86:       Vec            X,F
 87:       integer        ierr
 88:       type (userctx) user

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

 94:       write(*,*)"Inside FormFunction, user%xm=",user%xm

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

101:       call DAGetLocalVector(user%da,localX,ierr)
102:       call DAGlobalToLocalBegin(user%da,X,INSERT_VALUES,                &
103:      &     localX,ierr)
104:       call DAGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)

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

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

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

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

121: !  Restore vectors

123:       call VecRestoreArrayF90(localX,lx_v,ierr)
124:       call VecRestoreArrayF90(F,lf_v,ierr)

126: !  Insert values into global vector

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

131: !      call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
132: !      call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)

134:       return
135:       end subroutine formfunction
136:       end module f90module



140:       program main
141:       use f90module
142:       implicit none
143: !
144: !
145:  #include finclude/petscsys.h
146:  #include finclude/petscvec.h
147:  #include finclude/petscda.h
148:  #include finclude/petscis.h
149:  #include finclude/petscmat.h
150:  #include finclude/petscksp.h
151:  #include finclude/petscpc.h
152:  #include finclude/petscsnes.h
153: #include "finclude/petscvec.h90"
154: #include "finclude/petscda.h90"

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

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

180:       external FormInitialGuess,FormJacobian,FormLineSearch
181:       external FormPostCheck,FormPreCheck

183: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184: !  Initialize program
185: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

187:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
188:       call MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)

190: !  Initialize problem parameters

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


204: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205: !  Create nonlinear solver context
206: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

208:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

210: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: !  Create vector data structures; set function evaluation routine
212: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

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

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

236:       call DACreateGlobalVector(user%da,x,ierr)
237:       call VecDuplicate(x,r,ierr)

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

241:       call DAGetCorners(user%da,user%xs,user%ys,PETSC_NULL_INTEGER,     &
242:      &     user%xm,user%ym,PETSC_NULL_INTEGER,ierr)
243:       call DAGetGhostCorners(user%da,user%gxs,user%gys,                 &
244:      &     PETSC_NULL_INTEGER,user%gxm,user%gym,                        &
245:      &     PETSC_NULL_INTEGER,ierr)

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

250:       user%xs  = user%xs+1
251:       user%ys  = user%ys+1
252:       user%gxs = user%gxs+1
253:       user%gys = user%gys+1

255:       user%ye  = user%ys+user%ym-1
256:       user%xe  = user%xs+user%xm-1
257:       user%gye = user%gys+user%gym-1
258:       user%gxe = user%gxs+user%gxm-1

260: !  Set function evaluation routine and vector

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

264: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265: !  Create matrix data structure; set Jacobian evaluation routine
266: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

288:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_mf',         &
289:      &     matrix_free,ierr)
290:       if (matrix_free .eq. 0) then
291:         call DAGetMatrix(user%da,MATAIJ,J,ierr)
292:         call SNESSetJacobian(snes,J,J,FormJacobian,user,ierr)
293:       endif

295: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
296: !  Customize nonlinear solver; set runtime options
297: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

301:       call SNESSetFromOptions(snes,ierr)
302:       test_linesearch = 0
303:       test_check      = 0
304:       call PetscOptionsGetTruth(PETSC_NULL_CHARACTER,'-test_check',
305:      &                          test_check,PETSC_NULL_INTEGER,ierr)
306:       if (test_check.eq.1) then
307:         call SNESLineSearchSetPreCheck(snes,FormPreCheck,user,ierr)
308:         call SNESLineSearchSetPostCheck(snes,FormPostCheck,user,ierr)
309:       else
310:         call PetscOptionsGetTruth(PETSC_NULL_CHARACTER,
311:      &       '-test_linesearch',test_linesearch,PETSC_NULL_INTEGER,ierr)
312:         if (test_linesearch.eq.1) then
313:           call SNESLineSearchSet(snes,FormLineSearch,user,ierr)
314:         end if
315:       end if
316: 
317: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318: !  Evaluate initial guess; then solve nonlinear system.
319: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

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

336: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
337: !  Free work space.  All PETSc objects should be destroyed when they
338: !  are no longer needed.
339: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

341:       if (matrix_free .eq. 0) call MatDestroy(J,ierr)
342:       call VecDestroy(x,ierr)
343:       call VecDestroy(r,ierr)
344:       call SNESDestroy(snes,ierr)
345:       call DADestroy(user%da,ierr)
346:       call PetscFinalize(ierr)
347:       end

349: ! ---------------------------------------------------------------------
350: !
351: !  FormLineSearch - Applies the line search to the step size
352: !
353:       subroutine FormLineSearch(snes,user,x,f,g,y,w,fnorm,ynorm,gnorm,
354:      & flag,ierr)

356:       use f90module

358:  #include finclude/petscsys.h
359:  #include finclude/petscvec.h
360: #include "finclude/petscvec.h90"
361:  #include finclude/petscmat.h
362: #include "finclude/petscmat.h90"
363:  #include finclude/petscksp.h
364:  #include finclude/petscpc.h
365:  #include finclude/petscsnes.h

367:       SNES              snes
368:       type (userctx)    user
369:       Vec               x,f,g,y,w
370:       PetscReal fnorm,ynorm,gnorm
371:       PetscTruth           flag
372:       PetscErrorCode ierr

374:       PetscScalar       mone

376:         write(*,*)"Inside FormLineSearch, user%xm=",user%xm
377:       mone = -1.0d0
378:       flag = 0
379:       call VecNorm(y,NORM_2,ynorm,ierr)
380:       call VecAYPX(y,mone,x,ierr)
381:       call SNESComputeFunction(snes,y,g,ierr)
382:       call VecNorm(g,NORM_2,gnorm,ierr)
383:       return
384:       end

386: ! ---------------------------------------------------------------------
387: !
388: !  FormInitialGuess - Forms initial approximation.
389: !
390: !  Input Parameters:
391: !  X - vector
392: !
393: !  Output Parameter:
394: !  X - vector
395: !
396: !  Notes:
397: !  This routine serves as a wrapper for the lower-level routine
398: !  "InitialGuessLocal", where the actual computations are
399: !  done using the standard Fortran style of treating the local
400: !  vector data as a multidimensional array over the local mesh.
401: !  This routine merely handles ghost point scatters and accesses
402: !  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
403: !
404:       subroutine FormInitialGuess(user,X,ierr)
405:       use f90module
406:       implicit none

408: #include "finclude/petscvec.h90"
409:  #include finclude/petscsys.h
410:  #include finclude/petscvec.h
411:  #include finclude/petscda.h
412:  #include finclude/petscis.h
413:  #include finclude/petscmat.h
414:  #include finclude/petscksp.h
415:  #include finclude/petscpc.h
416:  #include finclude/petscsnes.h

418: !  Input/output variables:
419:       type (userctx)         user
420:       Vec      X
421:       integer  ierr
422: 
423: !  Declarations for use with local arrays:
424:       PetscScalar,pointer :: lx_v(:)
425:       Vec               localX

427:       0

429: !  Get a pointer to vector data.
430: !    - For default PETSc vectors, VecGetArray90() returns a pointer to
431: !      the data array. Otherwise, the routine is implementation dependent.
432: !    - You MUST call VecRestoreArrayF90() when you no longer need access to
433: !      the array.
434: !    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
435: !      and is useable from Fortran-90 Only.

437:       call DAGetLocalVector(user%da,localX,ierr)
438:       call VecGetArrayF90(localX,lx_v,ierr)

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

442:       call InitialGuessLocal(user,lx_v,ierr)

444: !  Restore vector

446:       call VecRestoreArrayF90(localX,lx_v,ierr)

448: !  Insert values into global vector

450:       call DALocalToGlobal(user%da,localX,INSERT_VALUES,X,ierr)
451:       call DARestoreLocalVector(user%da,localX,ierr)

453:       return
454:       end

456: ! ---------------------------------------------------------------------
457: !
458: !  InitialGuessLocal - Computes initial approximation, called by
459: !  the higher level routine FormInitialGuess().
460: !
461: !  Input Parameter:
462: !  x - local vector data
463: !
464: !  Output Parameters:
465: !  x - local vector data
466: !  ierr - error code
467: !
468: !  Notes:
469: !  This routine uses standard Fortran-style computations over a 2-dim array.
470: !
471:       subroutine InitialGuessLocal(user,x,ierr)
472:       use f90module
473:       implicit none

475:  #include finclude/petscsys.h
476:  #include finclude/petscvec.h
477:  #include finclude/petscda.h
478:  #include finclude/petscis.h
479:  #include finclude/petscmat.h
480:  #include finclude/petscksp.h
481:  #include finclude/petscpc.h
482:  #include finclude/petscsnes.h

484: !  Input/output variables:
485:       type (userctx)         user
486:       PetscScalar  x(user%gxs:user%gxe,                                         &
487:      &              user%gys:user%gye)
488:       integer ierr

490: !  Local variables:
491:       integer  i,j
492:       PetscScalar   temp1,temp,hx,hy
493:       PetscScalar   one

495: !  Set parameters

497:       0
498:       one    = 1.0
499:       hx     = one/(dble(user%mx-1))
500:       hy     = one/(dble(user%my-1))
501:       temp1  = user%lambda/(user%lambda + one)

503:       do 20 j=user%ys,user%ye
504:          temp = dble(min(j-1,user%my-j))*hy
505:          do 10 i=user%xs,user%xe
506:             if (i .eq. 1 .or. j .eq. 1                                  &
507:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
508:               x(i,j) = 0.0
509:             else
510:               x(i,j) = temp1 *                                          &
511:      &          sqrt(min(dble(min(i-1,user%mx-i)*hx),dble(temp)))
512:             endif
513:  10      continue
514:  20   continue

516:       return
517:       end

519: ! ---------------------------------------------------------------------
520: !
521: !  FormFunctionLocal - Computes nonlinear function, called by
522: !  the higher level routine FormFunction().
523: !
524: !  Input Parameter:
525: !  x - local vector data
526: !
527: !  Output Parameters:
528: !  f - local vector data, f(x)
529: !  ierr - error code
530: !
531: !  Notes:
532: !  This routine uses standard Fortran-style computations over a 2-dim array.
533: !
534:       subroutine FormFunctionLocal(x,f,user,ierr)
535:       use f90module

537:       implicit none

539: !  Input/output variables:
540:       type (userctx) user
541:       PetscScalar  x(user%gxs:user%gxe,                                         &
542:      &              user%gys:user%gye)
543:       PetscScalar  f(user%xs:user%xe,                                           &
544:      &              user%ys:user%ye)
545:       integer  ierr

547: !  Local variables:
548:       PetscScalar   two,one,hx,hy,hxdhy,hydhx,sc
549:       PetscScalar   u,uxx,uyy
550:       integer  i,j

552:       one    = 1.0
553:       two    = 2.0
554:       hx     = one/dble(user%mx-1)
555:       hy     = one/dble(user%my-1)
556:       sc     = hx*hy*user%lambda
557:       hxdhy  = hx/hy
558:       hydhx  = hy/hx

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

562:       do 20 j=user%ys,user%ye
563:          do 10 i=user%xs,user%xe
564:             if (i .eq. 1 .or. j .eq. 1                                  &
565:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
566:                f(i,j) = x(i,j)
567:             else
568:                u = x(i,j)
569:                uxx = hydhx * (two*u                                     &
570:      &                - x(i-1,j) - x(i+1,j))
571:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
572:                f(i,j) = uxx + uyy - sc*exp(u)
573:             endif
574:  10      continue
575:  20   continue

577:       return
578:       end

580: ! ---------------------------------------------------------------------
581: !
582: !  FormJacobian - Evaluates Jacobian matrix.
583: !
584: !  Input Parameters:
585: !  snes     - the SNES context
586: !  x        - input vector
587: !  dummy    - optional user-defined context, as set by SNESSetJacobian()
588: !             (not used here)
589: !
590: !  Output Parameters:
591: !  jac      - Jacobian matrix
592: !  jac_prec - optionally different preconditioning matrix (not used here)
593: !  flag     - flag indicating matrix structure
594: !
595: !  Notes:
596: !  This routine serves as a wrapper for the lower-level routine
597: !  "FormJacobianLocal", where the actual computations are
598: !  done using the standard Fortran style of treating the local
599: !  vector data as a multidimensional array over the local mesh.
600: !  This routine merely accesses the local vector data via
601: !  VecGetArrayF90() and VecRestoreArrayF90().
602: !
603: !  Notes:
604: !  Due to grid point reordering with DAs, we must always work
605: !  with the local grid points, and then transform them to the new
606: !  global numbering with the "ltog" mapping (via DAGetGlobalIndicesF90()).
607: !  We cannot work directly with the global numbers for the original
608: !  uniprocessor grid!
609: !
610: !  Two methods are available for imposing this transformation
611: !  when setting matrix entries:
612: !    (A) MatSetValuesLocal(), using the local ordering (including
613: !        ghost points!)
614: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
615: !        - Associate this map with the matrix by calling
616: !          MatSetLocalToGlobalMapping() once
617: !        - Set matrix entries using the local ordering
618: !          by calling MatSetValuesLocal()
619: !    (B) MatSetValues(), using the global ordering
620: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
621: !        - Then apply this map explicitly yourself
622: !        - Set matrix entries using the global ordering by calling
623: !          MatSetValues()
624: !  Option (A) seems cleaner/easier in many cases, and is the procedure
625: !  used in this example.
626: !
627:       subroutine FormJacobian(snes,X,jac,jac_prec,flag,user,ierr)
628:       use f90module
629:       implicit none

631:  #include finclude/petscsys.h
632:  #include finclude/petscvec.h
633:  #include finclude/petscda.h
634:  #include finclude/petscis.h
635:  #include finclude/petscmat.h
636:  #include finclude/petscksp.h
637:  #include finclude/petscpc.h
638:  #include finclude/petscsnes.h

640: #include "finclude/petscvec.h90"

642: !  Input/output variables:
643:       SNES         snes
644:       Vec          X
645:       Mat          jac,jac_prec
646:       MatStructure flag
647:       type(userctx) user
648:       integer      ierr

650: !  Declarations for use with local arrays:
651:       PetscScalar,pointer :: lx_v(:)
652:       Vec            localX

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

659:       call DAGetLocalVector(user%da,localX,ierr)
660:       call DAGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX,            &
661:      &     ierr)
662:       call DAGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)

664: !  Get a pointer to vector data

666:       call VecGetArrayF90(localX,lx_v,ierr)

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

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

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

677:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
678:       call VecRestoreArrayF90(localX,lx_v,ierr)
679:       call DARestoreLocalVector(user%da,localX,ierr)
680:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)

682: !  Set flag to indicate that the Jacobian matrix retains an identical
683: !  nonzero structure throughout all nonlinear iterations (although the
684: !  values of the entries change). Thus, we can save some work in setting
685: !  up the preconditioner (e.g., no need to redo symbolic factorization for
686: !  ILU/ICC preconditioners).
687: !   - If the nonzero structure of the matrix is different during
688: !     successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
689: !     must be used instead.  If you are unsure whether the matrix
690: !     structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
691: !   - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
692: !     believes your assertion and does not check the structure
693: !     of the matrix.  If you erroneously claim that the structure
694: !     is the same when it actually is not, the new preconditioner
695: !     will not function correctly.  Thus, use this optimization
696: !     feature with caution!

698:       flag = SAME_NONZERO_PATTERN

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

703: !       call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,ierr)

705:       return
706:       end

708: ! ---------------------------------------------------------------------
709: !
710: !  FormJacobianLocal - Computes Jacobian matrix, called by
711: !  the higher level routine FormJacobian().
712: !
713: !  Input Parameters:
714: !  x        - local vector data
715: !
716: !  Output Parameters:
717: !  jac      - Jacobian matrix
718: !  jac_prec - optionally different preconditioning matrix (not used here)
719: !  ierr     - error code
720: !
721: !  Notes:
722: !  This routine uses standard Fortran-style computations over a 2-dim array.
723: !
724: !  Notes:
725: !  Due to grid point reordering with DAs, we must always work
726: !  with the local grid points, and then transform them to the new
727: !  global numbering with the "ltog" mapping (via DAGetGlobalIndicesF90()).
728: !  We cannot work directly with the global numbers for the original
729: !  uniprocessor grid!
730: !
731: !  Two methods are available for imposing this transformation
732: !  when setting matrix entries:
733: !    (A) MatSetValuesLocal(), using the local ordering (including
734: !        ghost points!)
735: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
736: !        - Associate this map with the matrix by calling
737: !          MatSetLocalToGlobalMapping() once
738: !        - Set matrix entries using the local ordering
739: !          by calling MatSetValuesLocal()
740: !    (B) MatSetValues(), using the global ordering
741: !        - Use DAGetGlobalIndicesF90() to extract the local-to-global map
742: !        - Then apply this map explicitly yourself
743: !        - Set matrix entries using the global ordering by calling
744: !          MatSetValues()
745: !  Option (A) seems cleaner/easier in many cases, and is the procedure
746: !  used in this example.
747: !
748:       subroutine FormJacobianLocal(x,jac,jac_prec,user,ierr)
749:       use f90module
750:       implicit none

752:  #include finclude/petscsys.h
753:  #include finclude/petscvec.h
754:  #include finclude/petscda.h
755:  #include finclude/petscis.h
756:  #include finclude/petscmat.h
757:  #include finclude/petscksp.h
758:  #include finclude/petscpc.h
759:  #include finclude/petscsnes.h

761: !  Input/output variables:
762:       type (userctx) user
763:       PetscScalar  x(user%gxs:user%gxe,                                         &
764:      &              user%gys:user%gye)
765:       Mat      jac,jac_prec
766:       integer  ierr

768: !  Local variables:
769:       integer  row,col(5),i,j
770:       PetscScalar   two,one,hx,hy,hxdhy
771:       PetscScalar   hydhx,sc,v(5)

773: !  Set parameters

775:       one    = 1.0
776:       two    = 2.0
777:       hx     = one/dble(user%mx-1)
778:       hy     = one/dble(user%my-1)
779:       sc     = hx*hy
780:       hxdhy  = hx/hy
781:       hydhx  = hy/hx

783: !  Compute entries for the locally owned part of the Jacobian.
784: !   - Currently, all PETSc parallel matrix formats are partitioned by
785: !     contiguous chunks of rows across the processors.
786: !   - Each processor needs to insert only elements that it owns
787: !     locally (but any non-local elements will be sent to the
788: !     appropriate processor during matrix assembly).
789: !   - Here, we set all entries for a particular row at once.
790: !   - We can set matrix entries either using either
791: !     MatSetValuesLocal() or MatSetValues(), as discussed above.
792: !   - Note that MatSetValues() uses 0-based row and column numbers
793: !     in Fortran as well as in C.

795:       do 20 j=user%ys,user%ye
796:          row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
797:          do 10 i=user%xs,user%xe
798:             row = row + 1
799: !           boundary points
800:             if (i .eq. 1 .or. j .eq. 1                                  &
801:      &             .or. i .eq. user%mx .or. j .eq. user%my) then
802:                col(1) = row
803:                v(1)   = one
804:                call MatSetValuesLocal(jac,1,row,1,col,v,                &
805:      &                           INSERT_VALUES,ierr)
806: !           interior grid points
807:             else
808:                v(1) = -hxdhy
809:                v(2) = -hydhx
810:                v(3) = two*(hydhx + hxdhy)                               &
811:      &                  - sc*user%lambda*exp(x(i,j))
812:                v(4) = -hydhx
813:                v(5) = -hxdhy
814:                col(1) = row - user%gxm
815:                col(2) = row - 1
816:                col(3) = row
817:                col(4) = row + 1
818:                col(5) = row + user%gxm
819:                call MatSetValuesLocal(jac,1,row,5,col,v,                &
820:      &                                INSERT_VALUES,ierr)
821:             endif
822:  10      continue
823:  20   continue

825:       return
826:       end

828:       subroutine FormPreCheck(snes,X,Y,user,changed_Y,ierr)
829:       use f90module

831:       SNES           snes
832:       Vec            X,Y
833:       type (userctx) user
834:       PetscTruth     changed_Y
835:       PetscErrorCode ierr

837:              write(*,*)"Inside formPreCheck, user%xm=",user%xm
838:       end subroutine formPreCheck

840:       subroutine FormPostCheck(snes,X,Y,W,user,changed_Y,changed_W,ierr)
841:       use f90module

843:       SNES           snes
844:       Vec            X,Y,W
845:       type (userctx) user
846:       PetscTruth     changed_Y,changed_W
847:       PetscErrorCode ierr

849:              write(*,*)"Inside formPostCheck, user%xm=",user%xm
850:       end subroutine formPostCheck