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

 74:  #include include/finclude/petsc.h
 75:  #include include/finclude/petscvec.h
 76:  #include include/finclude/petscda.h
 77:  #include include/finclude/petscis.h
 78:  #include include/finclude/petscmat.h
 79:  #include include/finclude/petscksp.h
 80:  #include include/finclude/petscpc.h
 81:  #include include/finclude/petscsnes.h

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


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

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

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

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

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

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

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

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

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

123: !  Restore vectors

125:       call VecRestoreArrayF90(localX,lx_v,ierr)
126:       call VecRestoreArrayF90(F,lf_v,ierr)

128: !  Insert values into global vector

130:       call DARestoreLocalVector(user%da,localX,ierr)
131:       call PetscLogFlops(11*user%ym*user%xm,ierr)

133: !      call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
134: !      call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)

136:       return
137:       end subroutine formfunction
138:       end module f90module



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

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

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

182:       external FormInitialGuess,FormJacobian,FormLineSearch
183:       external FormPostCheck,FormPreCheck

185: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186: !  Initialize program
187: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

189:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
190:       call MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)

192: !  Initialize problem parameters

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


206: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207: !  Create nonlinear solver context
208: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

210:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

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

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

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

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

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

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

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

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

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

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

262: !  Set function evaluation routine and vector

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

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

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

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

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

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

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

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

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

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

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

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

358:       use f90module

360:  #include include/finclude/petsc.h
361:  #include include/finclude/petscvec.h
362: #include "include/finclude/petscvec.h90"
363:  #include include/finclude/petscmat.h
364: #include "include/finclude/petscmat.h90"
365:  #include include/finclude/petscksp.h
366:  #include include/finclude/petscpc.h
367:  #include include/finclude/petscsnes.h

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

376:       PetscScalar       mone

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

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

410: #include "include/finclude/petscvec.h90"
411:  #include include/finclude/petsc.h
412:  #include include/finclude/petscvec.h
413:  #include include/finclude/petscda.h
414:  #include include/finclude/petscis.h
415:  #include include/finclude/petscmat.h
416:  #include include/finclude/petscksp.h
417:  #include include/finclude/petscpc.h
418:  #include include/finclude/petscsnes.h

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

429:       0

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

439:       call DAGetLocalVector(user%da,localX,ierr)
440:       call VecGetArrayF90(localX,lx_v,ierr)

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

444:       call InitialGuessLocal(user,lx_v,ierr)

446: !  Restore vector

448:       call VecRestoreArrayF90(localX,lx_v,ierr)

450: !  Insert values into global vector

452:       call DALocalToGlobal(user%da,localX,INSERT_VALUES,X,ierr)
453:       call DARestoreLocalVector(user%da,localX,ierr)

455:       return
456:       end

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

477:  #include include/finclude/petsc.h
478:  #include include/finclude/petscvec.h
479:  #include include/finclude/petscda.h
480:  #include include/finclude/petscis.h
481:  #include include/finclude/petscmat.h
482:  #include include/finclude/petscksp.h
483:  #include include/finclude/petscpc.h
484:  #include include/finclude/petscsnes.h

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

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

497: !  Set parameters

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

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

518:       return
519:       end

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

539:       implicit none

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

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

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

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

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

579:       return
580:       end

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

633:  #include include/finclude/petsc.h
634:  #include include/finclude/petscvec.h
635:  #include include/finclude/petscda.h
636:  #include include/finclude/petscis.h
637:  #include include/finclude/petscmat.h
638:  #include include/finclude/petscksp.h
639:  #include include/finclude/petscpc.h
640:  #include include/finclude/petscsnes.h

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

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

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

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

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

666: !  Get a pointer to vector data

668:       call VecGetArrayF90(localX,lx_v,ierr)

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

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

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

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

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

700:       flag = SAME_NONZERO_PATTERN

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

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

707:       return
708:       end

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

754:  #include include/finclude/petsc.h
755:  #include include/finclude/petscvec.h
756:  #include include/finclude/petscda.h
757:  #include include/finclude/petscis.h
758:  #include include/finclude/petscmat.h
759:  #include include/finclude/petscksp.h
760:  #include include/finclude/petscpc.h
761:  #include include/finclude/petscsnes.h

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

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

775: !  Set parameters

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

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

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

827:       return
828:       end

830:       subroutine FormPreCheck(snes,X,Y,user,changed_Y,ierr)
831:       use f90module

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

839:              write(*,*)"Inside formPreCheck, user%xm=",user%xm
840:       end subroutine formPreCheck

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

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

851:              write(*,*)"Inside formPostCheck, user%xm=",user%xm
852:       end subroutine formPostCheck