Actual source code: ex5f.F

petsc-3.4.5 2014-06-29
  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:       implicit none
 35: !
 36: !  We place common blocks, variable declarations, and other include files
 37: !  needed for this code in the single file ex5f.h.  We then need to include
 38: !  only this file throughout the various routines in this program.  See
 39: !  additional comments in the file ex5f.h.
 40: !
 41: #include "ex5f.h"

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


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

 65:       external FormInitialGuess
 66:       external FormFunctionLocal,FormJacobianLocal

 68: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69: !  Initialize program
 70: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 72:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 73:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
 74:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 76: !  Initialize problem parameters

 78:       i1 = 1
 79:       i4 = -4
 80:       lambda_max = 6.81
 81:       lambda_min = 0.0
 82:       lambda     = 6.0
 83:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',lambda,                &
 84:      &                           flg,ierr)
 85:       if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
 86:          if (rank .eq. 0) write(6,*) 'Lambda is out of range'
 87:          SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
 88:       endif

 90: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 91: !  Create nonlinear solver context
 92: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 94:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

 96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97: !  Create vector data structures; set function evaluation routine
 98: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

102: ! This really needs only the star-type stencil, but we use the box
103: ! stencil temporarily.
104:       call DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,            &
105:      &     DMDA_BOUNDARY_NONE,                                          &
106:      &     DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE,i1,i1,                &
107:      &     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)

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

112:       call DMCreateGlobalVector(da,x,ierr)
113:       call VecDuplicate(x,r,ierr)

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

117:       call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,            &
118:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
119:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
120:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
121:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
122:      &               PETSC_NULL_INTEGER,ierr)
123:       call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,                      &
124:      &     PETSC_NULL_INTEGER,ierr)
125:       call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,             &
126:      &     PETSC_NULL_INTEGER,ierr)

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

131:       xs  = xs+1
132:       ys  = ys+1
133:       gxs = gxs+1
134:       gys = gys+1

136:       ye  = ys+ym-1
137:       xe  = xs+xm-1
138:       gye = gys+gym-1
139:       gxe = gxs+gxm-1

141: !  Set function evaluation routine and vector

143:       call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,            &
144:      &                              PETSC_NULL_OBJECT,ierr)
145:       call DMDASNESSetJacobianLocal(da,FormJacobianLocal,                          &
146:      &                              PETSC_NULL_OBJECT,ierr)
147:       call SNESSetDM(snes,da,ierr)

149: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: !  Customize nonlinear solver; set runtime options
151: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

155:           call SNESSetFromOptions(snes,ierr)
156: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: !  Evaluate initial guess; then solve nonlinear system.
158: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

165:       call FormInitialGuess(x,ierr)
166:       call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
167:       call SNESGetIterationNumber(snes,its,ierr)
168:       if (rank .eq. 0) then
169:          write(6,100) its
170:       endif
171:   100 format('Number of SNES iterations = ',i5)


174: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: !  Free work space.  All PETSc objects should be destroyed when they
176: !  are no longer needed.
177: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

179:       call VecDestroy(x,ierr)
180:       call VecDestroy(r,ierr)
181:       call SNESDestroy(snes,ierr)
182:       call DMDestroy(da,ierr)
183:       call PetscFinalize(ierr)
184:       end

186: ! ---------------------------------------------------------------------
187: !
188: !  FormInitialGuess - Forms initial approximation.
189: !
190: !  Input Parameters:
191: !  X - vector
192: !
193: !  Output Parameter:
194: !  X - vector
195: !
196: !  Notes:
197: !  This routine serves as a wrapper for the lower-level routine
198: !  "ApplicationInitialGuess", where the actual computations are
199: !  done using the standard Fortran style of treating the local
200: !  vector data as a multidimensional array over the local mesh.
201: !  This routine merely handles ghost point scatters and accesses
202: !  the local vector data via VecGetArray() and VecRestoreArray().
203: !
204:       subroutine FormInitialGuess(X,ierr)
205:       implicit none

207: #include "ex5f.h"

209: !  Input/output variables:
210:       Vec      X
211:       PetscErrorCode  ierr

213: !  Declarations for use with local arrays:
214:       PetscScalar lx_v(0:1)
215:       PetscOffset lx_i
216:       Vec         localX

218:       0

220: !  Get a pointer to vector data.
221: !    - For default PETSc vectors, VecGetArray() returns a pointer to
222: !      the data array.  Otherwise, the routine is implementation dependent.
223: !    - You MUST call VecRestoreArray() when you no longer need access to
224: !      the array.
225: !    - Note that the Fortran interface to VecGetArray() differs from the
226: !      C version.  See the users manual for details.

228:       call DMGetLocalVector(da,localX,ierr)
229:       call VecGetArray(localX,lx_v,lx_i,ierr)

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

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

235: !  Restore vector

237:       call VecRestoreArray(localX,lx_v,lx_i,ierr)

239: !  Insert values into global vector

241:       call DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X,ierr)
242:       call DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X,ierr)
243:       call DMRestoreLocalVector(da,localX,ierr)
244:       return
245:       end

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

265: #include "ex5f.h"

267: !  Input/output variables:
268:       PetscScalar  x(gxs:gxe,gys:gye)
269:       PetscErrorCode ierr

271: !  Local variables:
272:       PetscInt  i,j
273:       PetscReal temp1,temp,one,hx,hy

275: !  Set parameters

277:       0
278:       one    = 1.0
279:       hx     = one/((mx-1))
280:       hy     = one/((my-1))
281:       temp1  = lambda/(lambda + one)

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

296:       return
297:       end

299: ! ---------------------------------------------------------------------
300: !
301: !  FormFunctionLocal - Computes nonlinear function, called by
302: !  the higher level routine FormFunction().
303: !
304: !  Input Parameter:
305: !  x - local vector data
306: !
307: !  Output Parameters:
308: !  f - local vector data, f(x)
309: !  ierr - error code
310: !
311: !  Notes:
312: !  This routine uses standard Fortran-style computations over a 2-dim array.
313: !
314: !
315:       subroutine FormFunctionLocal(info,x,f,dummy,ierr)

317:       implicit none

319: #include "ex5f.h"

321: !  Input/output variables:
322:       DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)
323:       PetscScalar x(gxs:gxe,gys:gye)
324:       PetscScalar f(xs:xe,ys:ye)
325:       PetscErrorCode     ierr
326:       PetscObject dummy

328: !  Local variables:
329:       PetscScalar two,one,hx,hy
330:       PetscScalar hxdhy,hydhx,sc
331:       PetscScalar u,uxx,uyy
332:       PetscInt  i,j

334:       xs     = info(DMDA_LOCAL_INFO_XS)+1
335:       xe     = xs+info(DMDA_LOCAL_INFO_XM)-1
336:       ys     = info(DMDA_LOCAL_INFO_YS)+1
337:       ye     = ys+info(DMDA_LOCAL_INFO_YM)-1
338:       mx     = info(DMDA_LOCAL_INFO_MX)
339:       my     = info(DMDA_LOCAL_INFO_MY)

341:       one    = 1.0
342:       two    = 2.0
343:       hx     = one/(mx-1)
344:       hy     = one/(my-1)
345:       sc     = hx*hy*lambda
346:       hxdhy  = hx/hy
347:       hydhx  = hy/hx

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

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

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

368:       return
369:       end

371: ! ---------------------------------------------------------------------
372: !
373: !  FormJacobianLocal - Computes Jacobian matrix, called by
374: !  the higher level routine FormJacobian().
375: !
376: !  Input Parameters:
377: !  x        - local vector data
378: !
379: !  Output Parameters:
380: !  jac      - Jacobian matrix
381: !  jac_prec - optionally different preconditioning matrix (not used here)
382: !  ierr     - error code
383: !
384: !  Notes:
385: !  This routine uses standard Fortran-style computations over a 2-dim array.
386: !
387: !  Notes:
388: !  Due to grid point reordering with DMDAs, we must always work
389: !  with the local grid points, and then transform them to the new
390: !  global numbering with the "ltog" mapping (via DMDAGetGlobalIndices()).
391: !  We cannot work directly with the global numbers for the original
392: !  uniprocessor grid!
393: !
394: !  Two methods are available for imposing this transformation
395: !  when setting matrix entries:
396: !    (A) MatSetValuesLocal(), using the local ordering (including
397: !        ghost points!)
398: !        - Use DMDAGetGlobalIndices() to extract the local-to-global map
399: !        - Associate this map with the matrix by calling
400: !          MatSetLocalToGlobalMapping() once
401: !        - Set matrix entries using the local ordering
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,str,ierr)
412:       implicit none

414: #include "ex5f.h"

416: !  Input/output variables:
417:       PetscScalar x(gxs:gxe,gys:gye)
418:       Mat         A,jac
419:       MatStructure str
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:       str = SAME_NONZERO_PATTERN
491:       return
492:       end