Actual source code: ex5f.F90

petsc-3.9.0 2018-04-07
Report Typos and Errors
  1: !
  2: !  Description: This example 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 <param>, where <param> indicates the nonlinearity of the problem
  7: !       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
  8: !
  9: !
 10: !!/*T
 11: !  Concepts: SNES^parallel Bratu example
 12: !  Concepts: DMDA^using distributed arrays;
 13: !  Processors: n
 14: !T*/


 17: !
 18: !  --------------------------------------------------------------------------
 19: !
 20: !  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 21: !  the partial differential equation
 22: !
 23: !          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 24: !
 25: !  with boundary conditions
 26: !
 27: !           u = 0  for  x = 0, x = 1, y = 0, y = 1.
 28: !
 29: !  A finite difference approximation with the usual 5-point stencil
 30: !  is used to discretize the boundary value problem to obtain a nonlinear
 31: !  system of equations.
 32: !
 33: !  --------------------------------------------------------------------------

 35:       program main
 36:  #include <petsc/finclude/petscsnes.h>
 37:       use petscdmda
 38:       use petscsnes
 39:       implicit none
 40: !
 41: !  We place common blocks, variable declarations, and other include files
 42: !  needed for this code in the single file ex5f.h.  We then need to include
 43: !  only this file throughout the various routines in this program.  See
 44: !  additional comments in the file ex5f.h.
 45: !
 46: #include "ex5f.h"

 48: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49: !                   Variable declarations
 50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51: !
 52: !  Variables:
 53: !     snes        - nonlinear solver
 54: !     x, r        - solution, residual vectors
 55: !     its         - iterations for convergence
 56: !
 57: !  See additional variable declarations in the file ex5f.h
 58: !
 59:       SNES           snes
 60:       Vec            x,r
 61:       PetscInt       its,i1,i4
 62:       PetscErrorCode ierr
 63:       PetscReal      lambda_max,lambda_min
 64:       PetscBool      flg
 65:       DM             da

 67: !  Note: Any user-defined Fortran routines (such as FormJacobianLocal)
 68: !  MUST be declared as external.

 70:       external FormInitialGuess
 71:       external FormFunctionLocal,FormJacobianLocal
 72:       external MySNESConverged

 74: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75: !  Initialize program
 76: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 78:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 79:       if (ierr .ne. 0) then
 80:         print*,'Unable to initialize PETSc'
 81:         stop
 82:       endif
 83:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
 84:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 86: !  Initialize problem parameters

 88:       i1 = 1
 89:       i4 = 4
 90:       lambda_max = 6.81
 91:       lambda_min = 0.0
 92:       lambda     = 6.0
 93:       call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,flg,ierr)
 94:       if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then; SETERRA(PETSC_COMM_WORLD,1,'Lambda out of range'); endif

 96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97: !  Create nonlinear solver context
 98: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

100:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

102: !  Set convergence test routine if desired

104:       call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_snes_convergence',flg,ierr)
105:       if (flg) then
106:         call SNESSetConvergenceTest(snes,MySNESConverged,0,PETSC_NULL_FUNCTION,ierr)
107:       endif

109: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: !  Create vector data structures; set function evaluation routine
111: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

115: ! This really needs only the star-type stencil, but we use the box
116: ! stencil temporarily.
117:       call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,      &
118:      &     DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE,i1,i1,                &
119:      &     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
120:       call DMSetFromOptions(da,ierr)
121:       call DMSetUp(da,ierr)

123: !  Extract global and local vectors from DMDA; then duplicate for remaining
124: !  vectors that are the same types

126:       call DMCreateGlobalVector(da,x,ierr)
127:       call VecDuplicate(x,r,ierr)

129: !  Get local grid boundaries (for 2-dimensional DMDA)

131:       call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,            &
132:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
133:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
134:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
135:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
136:      &               PETSC_NULL_INTEGER,ierr)
137:       call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr)
138:       call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr)

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

143:       xs  = xs+1
144:       ys  = ys+1
145:       gxs = gxs+1
146:       gys = gys+1

148:       ye  = ys+ym-1
149:       xe  = xs+xm-1
150:       gye = gys+gym-1
151:       gxe = gxs+gxm-1

153: !  Set function evaluation routine and vector

155:       call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,da,ierr)
156:       call DMDASNESSetJacobianLocal(da,FormJacobianLocal,da,ierr)
157:       call SNESSetDM(snes,da,ierr)

159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: !  Customize nonlinear solver; set runtime options
161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

165:           call SNESSetFromOptions(snes,ierr)
166: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: !  Evaluate initial guess; then solve nonlinear system.
168: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

175:       call FormInitialGuess(x,ierr)
176:       call SNESSolve(snes,PETSC_NULL_VEC,x,ierr)
177:       call SNESGetIterationNumber(snes,its,ierr)
178:       if (rank .eq. 0) then
179:          write(6,100) its
180:       endif
181:   100 format('Number of SNES iterations = ',i5)


184: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: !  Free work space.  All PETSc objects should be destroyed when they
186: !  are no longer needed.
187: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

189:       call VecDestroy(x,ierr)
190:       call VecDestroy(r,ierr)
191:       call SNESDestroy(snes,ierr)
192:       call DMDestroy(da,ierr)
193:       call PetscFinalize(ierr)
194:       end

196: ! ---------------------------------------------------------------------
197: !
198: !  FormInitialGuess - Forms initial approximation.
199: !
200: !  Input Parameters:
201: !  X - vector
202: !
203: !  Output Parameter:
204: !  X - vector
205: !
206: !  Notes:
207: !  This routine serves as a wrapper for the lower-level routine
208: !  "ApplicationInitialGuess", where the actual computations are
209: !  done using the standard Fortran style of treating the local
210: !  vector data as a multidimensional array over the local mesh.
211: !  This routine merely handles ghost point scatters and accesses
212: !  the local vector data via VecGetArray() and VecRestoreArray().
213: !
214:       subroutine FormInitialGuess(X,ierr)
215:       use petscsnes
216:       implicit none

218: #include "ex5f.h"

220: !  Input/output variables:
221:       Vec      X
222:       PetscErrorCode  ierr

224: !  Declarations for use with local arrays:
225:       PetscScalar lx_v(0:1)
226:       PetscOffset lx_i

228:       0

230: !  Get a pointer to vector data.
231: !    - For default PETSc vectors, VecGetArray() returns a pointer to
232: !      the data array.  Otherwise, the routine is implementation dependent.
233: !    - You MUST call VecRestoreArray() when you no longer need access to
234: !      the array.
235: !    - Note that the Fortran interface to VecGetArray() differs from the
236: !      C version.  See the users manual for details.

238:       call VecGetArray(X,lx_v,lx_i,ierr)

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

242:       call InitialGuessLocal(lx_v(lx_i),ierr)

244: !  Restore vector

246:       call VecRestoreArray(X,lx_v,lx_i,ierr)

248:       return
249:       end

251: ! ---------------------------------------------------------------------
252: !
253: !  InitialGuessLocal - Computes initial approximation, called by
254: !  the higher level routine FormInitialGuess().
255: !
256: !  Input Parameter:
257: !  x - local vector data
258: !
259: !  Output Parameters:
260: !  x - local vector data
261: !  ierr - error code
262: !
263: !  Notes:
264: !  This routine uses standard Fortran-style computations over a 2-dim array.
265: !
266:       subroutine InitialGuessLocal(x,ierr)
267:       use petscsnes
268:       implicit none

270: #include "ex5f.h"

272: !  Input/output variables:
273:       PetscScalar    x(xs:xe,ys:ye)
274:       PetscErrorCode ierr

276: !  Local variables:
277:       PetscInt  i,j
278:       PetscReal temp1,temp,one,hx,hy

280: !  Set parameters

282:       0
283:       one    = 1.0
284:       hx     = one/((mx-1))
285:       hy     = one/((my-1))
286:       temp1  = lambda/(lambda + one)

288:       do 20 j=ys,ye
289:          temp = (min(j-1,my-j))*hy
290:          do 10 i=xs,xe
291:             if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
292:               x(i,j) = 0.0
293:             else
294:               x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,(temp)))
295:             endif
296:  10      continue
297:  20   continue

299:       return
300:       end

302: ! ---------------------------------------------------------------------
303: !
304: !  FormFunctionLocal - Computes nonlinear function, called by
305: !  the higher level routine FormFunction().
306: !
307: !  Input Parameter:
308: !  x - local vector data
309: !
310: !  Output Parameters:
311: !  f - local vector data, f(x)
312: !  ierr - error code
313: !
314: !  Notes:
315: !  This routine uses standard Fortran-style computations over a 2-dim array.
316: !
317: !
318:       subroutine FormFunctionLocal(info,x,f,da,ierr)
319:  #include <petsc/finclude/petscdmda.h>
320:       use petscsnes
321:       implicit none

323: #include "ex5f.h"
324:       DM da

326: !  Input/output variables:
327:       DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)
328:       PetscScalar x(gxs:gxe,gys:gye)
329:       PetscScalar f(xs:xe,ys:ye)
330:       PetscErrorCode     ierr

332: !  Local variables:
333:       PetscScalar two,one,hx,hy
334:       PetscScalar hxdhy,hydhx,sc
335:       PetscScalar u,uxx,uyy
336:       PetscInt  i,j

338:       xs     = info(DMDA_LOCAL_INFO_XS)+1
339:       xe     = xs+info(DMDA_LOCAL_INFO_XM)-1
340:       ys     = info(DMDA_LOCAL_INFO_YS)+1
341:       ye     = ys+info(DMDA_LOCAL_INFO_YM)-1
342:       mx     = info(DMDA_LOCAL_INFO_MX)
343:       my     = info(DMDA_LOCAL_INFO_MY)

345:       one    = 1.0
346:       two    = 2.0
347:       hx     = one/(mx-1)
348:       hy     = one/(my-1)
349:       sc     = hx*hy*lambda
350:       hxdhy  = hx/hy
351:       hydhx  = hy/hx

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

355:       do 20 j=ys,ye
356:          do 10 i=xs,xe
357:             if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
358:                f(i,j) = x(i,j)
359:             else
360:                u = x(i,j)
361:                uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j))
362:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
363:                f(i,j) = uxx + uyy - sc*exp(u)
364:             endif
365:  10      continue
366:  20   continue

368:       call PetscLogFlops(11.0d0*ym*xm,ierr)

370:       return
371:       end

373: ! ---------------------------------------------------------------------
374: !
375: !  FormJacobianLocal - Computes Jacobian matrix, called by
376: !  the higher level routine FormJacobian().
377: !
378: !  Input Parameters:
379: !  x        - local vector data
380: !
381: !  Output Parameters:
382: !  jac      - Jacobian matrix
383: !  jac_prec - optionally different preconditioning matrix (not used here)
384: !  ierr     - error code
385: !
386: !  Notes:
387: !  This routine uses standard Fortran-style computations over a 2-dim array.
388: !
389: !  Notes:
390: !  Due to grid point reordering with DMDAs, we must always work
391: !  with the local grid points, and then transform them to the new
392: !  global numbering with the "ltog" mapping
393: !  We cannot work directly with the global numbers for the original
394: !  uniprocessor grid!
395: !
396: !  Two methods are available for imposing this transformation
397: !  when setting matrix entries:
398: !    (A) MatSetValuesLocal(), using the local ordering (including
399: !        ghost points!)
400: !          by calling MatSetValuesLocal()
401: !    (B) MatSetValues(), using the global ordering
402: !        - Use DMDAGetGlobalIndices() to extract the local-to-global map
403: !        - Then apply this map explicitly yourself
404: !        - Set matrix entries using the global ordering by calling
405: !          MatSetValues()
406: !  Option (A) seems cleaner/easier in many cases, and is the procedure
407: !  used in this example.
408: !
409:       subroutine FormJacobianLocal(info,x,A,jac,da,ierr)
410:       use petscsnes
411:       implicit none

413: #include "ex5f.h"
414:       DM da
415: 
416: !  Input/output variables:
417:       PetscScalar x(gxs:gxe,gys:gye)
418:       Mat         A,jac
419:       PetscErrorCode  ierr
420:       DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)


423: !  Local variables:
424:       PetscInt  row,col(5),i,j,i1,i5
425:       PetscScalar two,one,hx,hy,v(5)
426:       PetscScalar hxdhy,hydhx,sc

428: !  Set parameters

430:       i1     = 1
431:       i5     = 5
432:       one    = 1.0
433:       two    = 2.0
434:       hx     = one/(mx-1)
435:       hy     = one/(my-1)
436:       sc     = hx*hy
437:       hxdhy  = hx/hy
438:       hydhx  = hy/hx

440: !  Compute entries for the locally owned part of the Jacobian.
441: !   - Currently, all PETSc parallel matrix formats are partitioned by
442: !     contiguous chunks of rows across the processors.
443: !   - Each processor needs to insert only elements that it owns
444: !     locally (but any non-local elements will be sent to the
445: !     appropriate processor during matrix assembly).
446: !   - Here, we set all entries for a particular row at once.
447: !   - We can set matrix entries either using either
448: !     MatSetValuesLocal() or MatSetValues(), as discussed above.
449: !   - Note that MatSetValues() uses 0-based row and column numbers
450: !     in Fortran as well as in C.

452:       do 20 j=ys,ye
453:          row = (j - gys)*gxm + xs - gxs - 1
454:          do 10 i=xs,xe
455:             row = row + 1
456: !           boundary points
457:             if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
458: !       Some f90 compilers need 4th arg to be of same type in both calls
459:                col(1) = row
460:                v(1)   = one
461:                call MatSetValuesLocal(jac,i1,row,i1,col,v,INSERT_VALUES,ierr)
462: !           interior grid points
463:             else
464:                v(1) = -hxdhy
465:                v(2) = -hydhx
466:                v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j))
467:                v(4) = -hydhx
468:                v(5) = -hxdhy
469:                col(1) = row - gxm
470:                col(2) = row - 1
471:                col(3) = row
472:                col(4) = row + 1
473:                col(5) = row + gxm
474:                call MatSetValuesLocal(jac,i1,row,i5,col,v, INSERT_VALUES,ierr)
475:             endif
476:  10      continue
477:  20   continue
478:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
479:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
480:       if (A .ne. jac) then
481:          call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
482:          call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
483:       endif
484:       return
485:       end

487: !
488: !     Simple convergence test based on the infinity norm of the residual being small
489: !
490:       subroutine MySNESConverged(snes,it,xnorm,snorm,fnorm,reason,dummy,ierr)
491:       use petscsnes
492:       implicit none

494:       SNES snes
495:       PetscInt it,dummy
496:       PetscReal xnorm,snorm,fnorm,nrm
497:       SNESConvergedReason reason
498:       Vec f
499:       PetscErrorCode ierr

501:       call SNESGetFunction(snes,f,PETSC_NULL_FUNCTION,dummy,ierr)
502:       call VecNorm(f,NORM_INFINITY,nrm,ierr)
503:       if (nrm .le. 1.e-5) reason = SNES_CONVERGED_FNORM_ABS

505:       end

507: !/*TEST
508: !
509: !   build:
510: !      requires: !complex !single
511: !
512: !   test:
513: !      nsize: 4
514: !      args: -snes_mf -da_processors_x 4 -da_processors_y 1 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
515: !
516: !   test:
517: !      suffix: 2
518: !      nsize: 4
519: !      args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
520: !
521: !   test:
522: !      suffix: 3
523: !      nsize: 3
524: !      args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
525: !
526: !   test:
527: !      suffix: 6
528: !      nsize: 1
529: !      args: -snes_monitor_short -my_snes_convergence
530: !
531: !TEST*/