Actual source code: ex5f.F

  1: ! "$Id: ex5f.F,v 1.80 2001/08/24 16:23:36 bsmith Exp $";
  2: !
  3: !  Description: This example solves a nonlinear system in parallel with SNES.
  4: !  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
  5: !  domain, using distributed arrays (DAs) to partition the parallel grid.
  6: !  The command line options include:
  7: !    -par <param>, where <param> indicates the nonlinearity of the problem
  8: !       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
  9: !
 10: !  Program usage:  mpirun -np <procs> ex5f [-help] [all PETSc options]
 11: !
 12: !/*T
 13: !  Concepts: SNES^parallel Bratu example
 14: !  Concepts: DA^using distributed arrays;
 15: !  Processors: n
 16: !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:       implicit none
 37: !
 38: !  We place common blocks, variable declarations, and other include files
 39: !  needed for this code in the single file ex5f.h.  We then need to include
 40: !  only this file throughout the various routines in this program.  See
 41: !  additional comments in the file ex5f.h.
 42: !
 43: #include "ex5f.h"

 45: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 46: !                   Variable declarations
 47: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 48: !
 49: !  Variables:
 50: !     snes        - nonlinear solver
 51: !     x, r        - solution, residual vectors
 52: !     J           - Jacobian matrix
 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:       Mat                    J,A
 60:       integer                its,flg,ierr
 61:       double precision       lambda_max,lambda_min
 62:       ISColoring             coloring
 63:       PetscTruth             adifor_jacobian,adiformf_jacobian

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

 68:       external FormInitialGuess
 69:       external FormFunctionLocal,FormJacobianLocal
 70: #if defined(PETSC_HAVE_ADIFOR)
 71:       external g_FormFunctionLocal,m_FormFunctionLocal
 72: #endif

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

 78:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 79:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
 80:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 82: !  Initialize problem parameters

 84:       lambda_max = 6.81
 85:       lambda_min = 0.0
 86:       lambda     = 6.0
 87:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',lambda,                &
 88:      &                           flg,ierr)
 89:       if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
 90:          if (rank .eq. 0) write(6,*) 'Lambda is out of range'
 91:          SETERRQ(1,' ',ierr)
 92:       endif

 94: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95: !  Create nonlinear solver context
 96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 98:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)

100: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: !  Create vector data structures; set function evaluation routine
102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

106: ! This really needs only the star-type stencil, but we use the box
107: ! stencil temporarily.
108:       call DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,-4,          &
109:      &     -4,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL_INTEGER,                   &
110:      &     PETSC_NULL_INTEGER,da,ierr)

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

115:       call DACreateGlobalVector(da,x,ierr)
116:       call VecDuplicate(x,r,ierr)

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

120:       call DAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,            &
121:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
122:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
123:      &               PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,                     &
124:      &               PETSC_NULL_INTEGER,ierr)
125:       call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,                      &
126:      &     PETSC_NULL_INTEGER,ierr)
127:       call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,             &
128:      &     PETSC_NULL_INTEGER,ierr)

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

133:       xs  = xs+1
134:       ys  = ys+1
135:       gxs = gxs+1
136:       gys = gys+1

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

143: !  Set function evaluation routine and vector

145:       call DASetLocalFunction(da,FormFunctionLocal,ierr)
146:       call DASetLocalJacobian(da,FormJacobianLocal,ierr)
147: #if defined(PETSC_HAVE_ADIFOR)
148:       call DASetLocalAdiforFunction(da,                                             &
149:      &         g_FormFunctionLocal,ierr)
150: #endif
151:       call SNESSetFunction(snes,r,SNESDAFormFunction,da,ierr)

153: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: !  Create matrix data structure; set Jacobian evaluation routine
155: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

157: !  Set Jacobian matrix data structure and default Jacobian evaluation
158: !  routine. User can override with:
159: !     -snes_mf : matrix-free Newton-Krylov method with no preconditioning
160: !                (unless user explicitly sets preconditioner)
161: !     -snes_mf_operator : form preconditioning matrix as set by the user,
162: !                         but use matrix-free approx for Jacobian-vector
163: !                         products within Newton-Krylov method
164: !

166:       call DAGetMatrix(da,MATMPIAIJ,J,ierr)

168: #if defined(PETSC_HAVE_ADIFOR)
169:       call PetscOptionsGetLogical(PETSC_NULL_CHARACTER                            &
170:      &     ,'-adiformf_jacobian',                                                 &
171:      &     adiformf_jacobian,PETSC_NULL_INTEGER,ierr)
172:       if (adiformf_jacobian .eq. 1) then
173:         call DASetLocalAdiforMFFunction(da,                                       &
174:      &         m_FormFunctionLocal,ierr)
175:         call MatRegisterDAAD(ierr)
176:         call MatCreateDAAD(da,A,ierr)
177:         call MatDAADSetSNES(A,snes,ierr)
178:       else
179:         A = J
180:       endif
181: #else 
182:       A = J
183: #endif

185:       call SNESSetJacobian(snes,A,J,SNESDAComputeJacobian,                        &
186:      &                       da,ierr)

188: #if defined(PETSC_HAVE_ADIFOR)
189:       call PetscOptionsGetLogical(PETSC_NULL_CHARACTER                            &
190:      &     ,'-adifor_jacobian',                                                   &
191:      &     adifor_jacobian,PETSC_NULL_INTEGER,ierr)
192:       if (adifor_jacobian .eq. 1) then
193:         call DAGetColoring(da,IS_COLORING_GHOSTED,                                 &
194:      &                     coloring,ierr)
195:         call MatSetColoring(J,coloring,ierr)
196:         call SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdifor,            &
197:      &                         da,ierr)
198:         call ISColoringDestroy(coloring,ierr)
199:       endif
200: #endif


203: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204: !  Customize nonlinear solver; set runtime options
205: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

209:           call SNESSetFromOptions(snes,ierr)

211: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212: !  Evaluate initial guess; then solve nonlinear system.
213: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

220:       call FormInitialGuess(x,ierr)
221:       call SNESSolve(snes,x,its,ierr)
222:       if (rank .eq. 0) then
223:          write(6,100) its
224:       endif
225:   100 format('Number of Newton iterations = ',i5)


228: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229: !  Free work space.  All PETSc objects should be destroyed when they
230: !  are no longer needed.
231: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

233:       if (A .ne. J) call MatDestroy(A,ierr)
234:       call MatDestroy(J,ierr)
235:       call VecDestroy(x,ierr)
236:       call VecDestroy(r,ierr)
237:       call SNESDestroy(snes,ierr)
238:       call DADestroy(da,ierr)
239:       call PetscFinalize(ierr)
240:       end

242: ! ---------------------------------------------------------------------
243: !
244: !  FormInitialGuess - Forms initial approximation.
245: !
246: !  Input Parameters:
247: !  X - vector
248: !
249: !  Output Parameter:
250: !  X - vector
251: !
252: !  Notes:
253: !  This routine serves as a wrapper for the lower-level routine
254: !  "ApplicationInitialGuess", where the actual computations are
255: !  done using the standard Fortran style of treating the local
256: !  vector data as a multidimensional array over the local mesh.
257: !  This routine merely handles ghost point scatters and accesses
258: !  the local vector data via VecGetArray() and VecRestoreArray().
259: !
260:       subroutine FormInitialGuess(X,ierr)
261:       implicit none

263: #include "ex5f.h"

265: !  Input/output variables:
266:       Vec      X
267:       integer  ierr

269: !  Declarations for use with local arrays:
270:       PetscScalar lx_v(0:1)
271:       PetscOffset lx_i
272:       Vec         localX

274:       0

276: !  Get a pointer to vector data.
277: !    - For default PETSc vectors, VecGetArray() returns a pointer to
278: !      the data array.  Otherwise, the routine is implementation dependent.
279: !    - You MUST call VecRestoreArray() when you no longer need access to
280: !      the array.
281: !    - Note that the Fortran interface to VecGetArray() differs from the
282: !      C version.  See the users manual for details.

284:       call DAGetLocalVector(da,localX,ierr)
285:       call VecGetArray(localX,lx_v,lx_i,ierr)

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

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

291: !  Restore vector

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

295: !  Insert values into global vector

297:       call DALocalToGlobal(da,localX,INSERT_VALUES,X,ierr)
298:       call DARestoreLocalVector(da,localX,ierr)
299:       return
300:       end

302: ! ---------------------------------------------------------------------
303: !
304: !  InitialGuessLocal - Computes initial approximation, called by
305: !  the higher level routine FormInitialGuess().
306: !
307: !  Input Parameter:
308: !  x - local vector data
309: !
310: !  Output Parameters:
311: !  x - local vector data
312: !  ierr - error code
313: !
314: !  Notes:
315: !  This routine uses standard Fortran-style computations over a 2-dim array.
316: !
317:       subroutine InitialGuessLocal(x,ierr)
318:       implicit none

320: #include "ex5f.h"

322: !  Input/output variables:
323:       PetscScalar  x(gxs:gxe,gys:gye)
324:       integer ierr

326: !  Local variables:
327:       integer  i,j,hxdhy,hydhx
328:       PetscScalar   temp1,temp,hx,hy,sc,one

330: !  Set parameters

332:       ierr   = 0
333:       one    = 1.0
334:       hx     = one/(dble(mx-1))
335:       hy     = one/(dble(my-1))
336:       sc     = hx*hy*lambda
337:       hxdhy  = hx/hy
338:       hydhx  = hy/hx
339:       temp1  = lambda/(lambda + one)

341:       do 20 j=ys,ye
342:          temp = dble(min(j-1,my-j))*hy
343:          do 10 i=xs,xe
344:             if (i .eq. 1 .or. j .eq. 1                                  &
345:      &             .or. i .eq. mx .or. j .eq. my) then
346:               x(i,j) = 0.0
347:             else
348:               x(i,j) = temp1 *                                          &
349:      &          sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
350:             endif
351:  10      continue
352:  20   continue

354:       return
355:       end

357: ! ---------------------------------------------------------------------
358: !
359: !  FormFunctionLocal - Computes nonlinear function, called by
360: !  the higher level routine FormFunction().
361: !
362: !  Input Parameter:
363: !  x - local vector data
364: !
365: !  Output Parameters:
366: !  f - local vector data, f(x)
367: !  ierr - error code
368: !
369: !  Notes:
370: !  This routine uses standard Fortran-style computations over a 2-dim array.
371: !
372: !     Process adifor: FormFunctionLocal
373: !
374:       subroutine FormFunctionLocal(info,x,f,dummy,ierr)

376:       implicit none

378: #include "ex5f.h"

380: !  Input/output variables:
381:       DALocalInfo info(DA_LOCAL_INFO_SIZE)
382:       PetscScalar x(gxs:gxe,gys:gye)
383:       PetscScalar f(xs:xe,ys:ye)
384:       integer     ierr
385:       PetscObject dummy

387: !  Local variables:
388:       PetscScalar   two,one,hx,hy,hxdhy,hydhx,sc
389:       PetscScalar   u,uxx,uyy
390:       integer  i,j


393:       xs     = info(DA_LOCAL_INFO_XS)+1
394:       xe     = xs+info(DA_LOCAL_INFO_XM)-1
395:       ys     = info(DA_LOCAL_INFO_YS)+1
396:       ye     = ys+info(DA_LOCAL_INFO_YM)-1
397:       mx     = info(DA_LOCAL_INFO_MX)
398:       my     = info(DA_LOCAL_INFO_MY)

400:       one    = 1.0
401:       two    = 2.0
402:       hx     = one/dble(mx-1)
403:       hy     = one/dble(my-1)
404:       sc     = hx*hy*lambda
405:       hxdhy  = hx/hy
406:       hydhx  = hy/hx

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

410:       do 20 j=ys,ye
411:          do 10 i=xs,xe
412:             if (i .eq. 1 .or. j .eq. 1                                  &
413:      &             .or. i .eq. mx .or. j .eq. my) then
414:                f(i,j) = x(i,j)
415:             else
416:                u = x(i,j)
417:                uxx = hydhx * (two*u                                     &
418:      &                - x(i-1,j) - x(i+1,j))
419:                uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
420:                f(i,j) = uxx + uyy - sc*exp(u)
421:             endif
422:  10      continue
423:  20   continue

425:       call PetscLogFlops(11*ym*xm,ierr)

427:       return
428:       end

430: ! ---------------------------------------------------------------------
431: !
432: !  FormJacobianLocal - Computes Jacobian matrix, called by
433: !  the higher level routine FormJacobian().
434: !
435: !  Input Parameters:
436: !  x        - local vector data
437: !
438: !  Output Parameters:
439: !  jac      - Jacobian matrix
440: !  jac_prec - optionally different preconditioning matrix (not used here)
441: !  ierr     - error code
442: !
443: !  Notes:
444: !  This routine uses standard Fortran-style computations over a 2-dim array.
445: !
446: !  Notes:
447: !  Due to grid point reordering with DAs, we must always work
448: !  with the local grid points, and then transform them to the new
449: !  global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
450: !  We cannot work directly with the global numbers for the original
451: !  uniprocessor grid!
452: !
453: !  Two methods are available for imposing this transformation
454: !  when setting matrix entries:
455: !    (A) MatSetValuesLocal(), using the local ordering (including
456: !        ghost points!)
457: !        - Use DAGetGlobalIndices() to extract the local-to-global map
458: !        - Associate this map with the matrix by calling
459: !          MatSetLocalToGlobalMapping() once
460: !        - Set matrix entries using the local ordering
461: !          by calling MatSetValuesLocal()
462: !    (B) MatSetValues(), using the global ordering
463: !        - Use DAGetGlobalIndices() to extract the local-to-global map
464: !        - Then apply this map explicitly yourself
465: !        - Set matrix entries using the global ordering by calling
466: !          MatSetValues()
467: !  Option (A) seems cleaner/easier in many cases, and is the procedure
468: !  used in this example.
469: !
470:       subroutine FormJacobianLocal(info,x,jac,ctx,ierr)
471:       implicit none

473: #include "ex5f.h"

475: !  Input/output variables:
476:       PetscScalar x(gxs:gxe,gys:gye)
477:       Mat         jac
478:       integer     ierr,ctx
479:       DALocalInfo info(DA_LOCAL_INFO_SIZE)
480: 

482: !  Local variables:
483:       integer  row,col(5),i,j
484:       PetscScalar   two,one,hx,hy,hxdhy,hydhx,sc,v(5)

486: !  Set parameters

488:       one    = 1.0
489:       two    = 2.0
490:       hx     = one/dble(mx-1)
491:       hy     = one/dble(my-1)
492:       sc     = hx*hy
493:       hxdhy  = hx/hy
494:       hydhx  = hy/hx

496: !  Compute entries for the locally owned part of the Jacobian.
497: !   - Currently, all PETSc parallel matrix formats are partitioned by
498: !     contiguous chunks of rows across the processors.
499: !   - Each processor needs to insert only elements that it owns
500: !     locally (but any non-local elements will be sent to the
501: !     appropriate processor during matrix assembly).
502: !   - Here, we set all entries for a particular row at once.
503: !   - We can set matrix entries either using either
504: !     MatSetValuesLocal() or MatSetValues(), as discussed above.
505: !   - Note that MatSetValues() uses 0-based row and column numbers
506: !     in Fortran as well as in C.

508:       do 20 j=ys,ye
509:          row = (j - gys)*gxm + xs - gxs - 1
510:          do 10 i=xs,xe
511:             row = row + 1
512: !           boundary points
513:             if (i .eq. 1 .or. j .eq. 1                                  &
514:      &             .or. i .eq. mx .or. j .eq. my) then
515: !       Some f90 compilers need 4th arg to be of same type in both calls
516:                col(1) = row
517:                v(1)   = one
518:                call MatSetValuesLocal(jac,1,row,1,col,v,                &
519:      &                           INSERT_VALUES,ierr)
520: !           interior grid points
521:             else
522:                v(1) = -hxdhy
523:                v(2) = -hydhx
524:                v(3) = two*(hydhx + hxdhy)                               &
525:      &                  - sc*lambda*exp(x(i,j))
526:                v(4) = -hydhx
527:                v(5) = -hxdhy
528:                col(1) = row - gxm
529:                col(2) = row - 1
530:                col(3) = row
531:                col(4) = row + 1
532:                col(5) = row + gxm
533:                call MatSetValuesLocal(jac,1,row,5,col,v,                &
534:      &                                INSERT_VALUES,ierr)
535:             endif
536:  10      continue
537:  20   continue
538:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
539:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
540:       return
541:       end