Actual source code: ex5f.F90

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

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

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


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

 68:       external FormInitialGuess
 69:       external FormFunctionLocal,FormJacobianLocal

 71: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 72: !  Initialize program
 73: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

 83: !  Initialize problem parameters

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

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

100:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: !  Create vector data structures; set function evaluation routine
104: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

108: ! This really needs only the star-type stencil, but we use the box
109: ! stencil temporarily.
110:       call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,              &
111:      &     DM_BOUNDARY_NONE,                                            &
112:      &     DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE,i1,i1,                &
113:      &     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
114:       call DMSetFromOptions(da,ierr)
115:       call DMSetUp(da,ierr)

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

120:       call DMCreateGlobalVector(da,x,ierr)
121:       call VecDuplicate(x,r,ierr)

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

125:       call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,            &
126:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
127:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
128:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
129:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
130:      &               PETSC_NULL_INTEGER,ierr)
131:       call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,                      &
132:      &     PETSC_NULL_INTEGER,ierr)
133:       call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,             &
134:      &     PETSC_NULL_INTEGER,ierr)

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

139:       xs  = xs+1
140:       ys  = ys+1
141:       gxs = gxs+1
142:       gys = gys+1

144:       ye  = ys+ym-1
145:       xe  = xs+xm-1
146:       gye = gys+gym-1
147:       gxe = gxs+gxm-1

149: !  Set function evaluation routine and vector

151:       call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,            &
152:      &                              0,ierr)
153:       call DMDASNESSetJacobianLocal(da,FormJacobianLocal,                          &
154:      &                              0,ierr)
155:       call SNESSetDM(snes,da,ierr)

157: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: !  Customize nonlinear solver; set runtime options
159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

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

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

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


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

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

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

216: #include "ex5f.h"

218: !  Input/output variables:
219:       Vec      X
220:       PetscErrorCode  ierr

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

226:       0

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

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

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

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

242: !  Restore vector

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

246:       return
247:       end

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

268: #include "ex5f.h"

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

274: !  Local variables:
275:       PetscInt  i,j
276:       PetscReal temp1,temp,one,hx,hy

278: !  Set parameters

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

286:       do 20 j=ys,ye
287:          temp = (min(j-1,my-j))*hy
288:          do 10 i=xs,xe
289:             if (i .eq. 1 .or. j .eq. 1                                  &
290:      &             .or. i .eq. mx .or. j .eq. my) then
291:               x(i,j) = 0.0
292:             else
293:               x(i,j) = temp1 *                                          &
294:      &          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,dummy,ierr)
319:  #include <petsc/finclude/petscdmda.h>
320:       use petscsnes
321:       implicit none

323: #include "ex5f.h"

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

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                                  &
358:      &             .or. i .eq. mx .or. j .eq. my) then
359:                f(i,j) = x(i,j)
360:             else
361:                u = x(i,j)
362:                uxx = hydhx * (two*u                                     &
363:      &                - x(i-1,j) - x(i+1,j))
364:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
365:                f(i,j) = uxx + uyy - sc*exp(u)
366:             endif
367:  10      continue
368:  20   continue

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

372:       return
373:       end

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

415: #include "ex5f.h"

417: !  Input/output variables:
418:       PetscScalar x(gxs:gxe,gys:gye)
419:       Mat         A,jac
420:       PetscErrorCode  ierr
421:       integer ctx
422:       DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)


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

430: !  Set parameters

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

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

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