Actual source code: ex14f.F90

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: !
  2: !
  3: !  Solves a nonlinear system in parallel with a user-defined
  4: !  Newton method that uses KSP to solve the linearized Newton sytems.  This solver
  5: !  is a very simplistic inexact Newton method.  The intent of this code is to
  6: !  demonstrate the repeated solution of linear sytems with the same nonzero pattern.
  7: !
  8: !  This is NOT the recommended approach for solving nonlinear problems with PETSc!
  9: !  We urge users to employ the SNES component for solving nonlinear problems whenever
 10: !  possible, as it offers many advantages over coding nonlinear solvers independently.
 11: !
 12: !  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
 13: !  domain, using distributed arrays (DMDAs) to partition the parallel grid.
 14: !
 15: !  The command line options include:
 16: !  -par <parameter>, where <parameter> indicates the problem's nonlinearity
 17: !     problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
 18: !  -mx <xg>, where <xg> = number of grid points in the x-direction
 19: !  -my <yg>, where <yg> = number of grid points in the y-direction
 20: !  -Nx <npx>, where <npx> = number of processors in the x-direction
 21: !  -Ny <npy>, where <npy> = number of processors in the y-direction
 22: !  -mf use matrix free for matrix vector product
 23: !
 24: !!/*T
 25: !   Concepts: KSP^writing a user-defined nonlinear solver
 26: !   Concepts: DMDA^using distributed arrays
 27: !   Processors: n
 28: !T*/


 31: !  ------------------------------------------------------------------------
 32: !
 33: !    Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 34: !    the partial differential equation
 35: !
 36: !            -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 37: !
 38: !    with boundary conditions
 39: !
 40: !             u = 0  for  x = 0, x = 1, y = 0, y = 1.
 41: !
 42: !    A finite difference approximation with the usual 5-point stencil
 43: !    is used to discretize the boundary value problem to obtain a nonlinear
 44: !    system of equations.
 45: !
 46: !    The SNES version of this problem is:  snes/examples/tutorials/ex5f.F
 47: !
 48: !  -------------------------------------------------------------------------
 49:       module mymoduleex14f
 50:  #include <petsc/finclude/petscksp.h>
 51:       use petscdmda
 52:       use petscksp
 53:       Vec      localX
 54:       PetscInt mx,my
 55:       Mat B
 56:       DM da
 57:       end module

 59:       program main
 60:       use mymoduleex14f
 61:       implicit none

 63:       MPI_Comm comm
 64:       Vec      X,Y,F
 65:       Mat      J
 66:       KSP      ksp

 68:       PetscInt  Nx,Ny,N,ifive,ithree
 69:       PetscBool  flg,nooutput,usemf
 70: !
 71: !      This is the routine to use for matrix-free approach
 72: !
 73:       external mymult

 75: !     --------------- Data to define nonlinear solver --------------
 76:       PetscReal   rtol,ttol
 77:       PetscReal   fnorm,ynorm,xnorm
 78:       PetscInt            max_nonlin_its,one
 79:       PetscInt            lin_its
 80:       PetscInt           i,m
 81:       PetscScalar        mone
 82:       PetscErrorCode ierr

 84:       mone           = -1.0
 85:       rtol           = 1.e-8
 86:       max_nonlin_its = 10
 87:       one            = 1
 88:       ifive          = 5
 89:       ithree         = 3

 91:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 92:       if (ierr .ne. 0) then
 93:         print*,'Unable to initialize PETSc'
 94:         stop
 95:       endif
 96:       comm = PETSC_COMM_WORLD

 98: !  Initialize problem parameters

100: !
101:       mx = 4
102:       my = 4
103:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr)
104:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my',my,flg,ierr)
105:       N = mx*my

107:       nooutput = .false.
108:       call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-no_output',nooutput,ierr)

110: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: !     Create linear solver context
112: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

114:       call KSPCreate(comm,ksp,ierr)

116: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: !     Create vector data structures
118: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

120: !
121: !  Create distributed array (DMDA) to manage parallel grid and vectors
122: !
123:       Nx = PETSC_DECIDE
124:       Ny = PETSC_DECIDE
125:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-Nx',Nx,flg,ierr)
126:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-Ny',Ny,flg,ierr)
127:       call DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,mx,my,Nx,Ny,one,one,          &
128:      &                  PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
129:       call DMSetFromOptions(da,ierr)
130:       call DMSetUp(da,ierr)
131: !
132: !  Extract global and local vectors from DMDA then duplicate for remaining
133: !  vectors that are the same types
134: !
135:        call DMCreateGlobalVector(da,X,ierr)
136:        call DMCreateLocalVector(da,localX,ierr)
137:        call VecDuplicate(X,F,ierr)
138:        call VecDuplicate(X,Y,ierr)


141: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: !     Create matrix data structure for Jacobian
143: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: !
145: !     Note:  For the parallel case, vectors and matrices MUST be partitioned
146: !     accordingly.  When using distributed arrays (DMDAs) to create vectors,
147: !     the DMDAs determine the problem partitioning.  We must explicitly
148: !     specify the local matrix dimensions upon its creation for compatibility
149: !     with the vector distribution.
150: !
151: !     Note: Here we only approximately preallocate storage space for the
152: !     Jacobian.  See the users manual for a discussion of better techniques
153: !     for preallocating matrix memory.
154: !
155:       call VecGetLocalSize(X,m,ierr)
156:       call MatCreateAIJ(comm,m,m,N,N,ifive,PETSC_NULL_INTEGER,ithree,PETSC_NULL_INTEGER,B,ierr)

158: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: !     if usemf is on then matrix vector product is done via matrix free
160: !     approach. Note this is just an example, and not realistic because
161: !     we still use the actual formed matrix, but in reality one would
162: !     provide their own subroutine that would directly do the matrix
163: !     vector product and not call MatMult()
164: !     Note: we put B into a module so it will be visible to the
165: !     mymult() routine
166:       usemf = .false.
167:       call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-mf',usemf,ierr)
168:       if (usemf) then
169:          call MatCreateShell(comm,m,m,N,N,PETSC_NULL_INTEGER,J,ierr)
170:          call MatShellSetOperation(J,MATOP_MULT,mymult,ierr)
171:       else
172: !        If not doing matrix free then matrix operator, J,  and matrix used
173: !        to construct preconditioner, B, are the same
174:         J = B
175:       endif

177: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: !     Customize linear solver set runtime options
179: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180: !
181: !     Set runtime options (e.g., -ksp_monitor -ksp_rtol <rtol> -ksp_type <type>)
182: !
183:        call KSPSetFromOptions(ksp,ierr)

185: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186: !     Evaluate initial guess
187: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

189:        call FormInitialGuess(X,ierr)
190:        call ComputeFunction(X,F,ierr)
191:        call VecNorm(F,NORM_2,fnorm,ierr)
192:        ttol = fnorm*rtol
193:        if (.not. nooutput) then
194:          print*, 'Initial function norm ',fnorm
195:        endif

197: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: !     Solve nonlinear system with a user-defined method
199: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

201: !  This solver is a very simplistic inexact Newton method, with no
202: !  no damping strategies or bells and whistles. The intent of this code
203: !  is merely to demonstrate the repeated solution with KSP of linear
204: !  sytems with the same nonzero structure.
205: !
206: !  This is NOT the recommended approach for solving nonlinear problems
207: !  with PETSc!  We urge users to employ the SNES component for solving
208: !  nonlinear problems whenever possible with application codes, as it
209: !  offers many advantages over coding nonlinear solvers independently.

211:        do 10 i=0,max_nonlin_its

213: !  Compute the Jacobian matrix.  See the comments in this routine for
214: !  important information about setting the flag mat_flag.

216:          call ComputeJacobian(X,B,ierr)

218: !  Solve J Y = F, where J is the Jacobian matrix.
219: !    - First, set the KSP linear operators.  Here the matrix that
220: !      defines the linear system also serves as the preconditioning
221: !      matrix.
222: !    - Then solve the Newton system.

224:          call KSPSetOperators(ksp,J,B,ierr)
225:          call KSPSolve(ksp,F,Y,ierr)

227: !  Compute updated iterate

229:          call VecNorm(Y,NORM_2,ynorm,ierr)
230:          call VecAYPX(Y,mone,X,ierr)
231:          call VecCopy(Y,X,ierr)
232:          call VecNorm(X,NORM_2,xnorm,ierr)
233:          call KSPGetIterationNumber(ksp,lin_its,ierr)
234:          if (.not. nooutput) then
235:            print*,'linear solve iterations = ',lin_its,' xnorm = ',xnorm,' ynorm = ',ynorm
236:          endif

238: !  Evaluate nonlinear function at new location

240:          call ComputeFunction(X,F,ierr)
241:          call VecNorm(F,NORM_2,fnorm,ierr)
242:          if (.not. nooutput) then
243:            print*, 'Iteration ',i+1,' function norm',fnorm
244:          endif

246: !  Test for convergence

248:        if (fnorm .le. ttol) then
249:          if (.not. nooutput) then
250:            print*,'Converged: function norm ',fnorm,' tolerance ',ttol
251:          endif
252:          goto 20
253:        endif
254:  10   continue
255:  20   continue

257:       write(6,100) i+1
258:  100  format('Number of SNES iterations =',I2)

260: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
261: !     Free work space.  All PETSc objects should be destroyed when they
262: !     are no longer needed.
263: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

265:        call MatDestroy(B,ierr)
266:        if (usemf) then
267:          call MatDestroy(J,ierr)
268:        endif
269:        call VecDestroy(localX,ierr)
270:        call VecDestroy(X,ierr)
271:        call VecDestroy(Y,ierr)
272:        call VecDestroy(F,ierr)
273:        call KSPDestroy(ksp,ierr)
274:        call DMDestroy(da,ierr)
275:        call PetscFinalize(ierr)
276:        end

278: ! -------------------------------------------------------------------
279: !
280: !   FormInitialGuess - Forms initial approximation.
281: !
282: !   Input Parameters:
283: !   X - vector
284: !
285: !   Output Parameter:
286: !   X - vector
287: !
288:       subroutine FormInitialGuess(X,ierr)
289:       use mymoduleex14f
290:       implicit none

292:       PetscErrorCode    ierr
293:       PetscOffset      idx
294:       Vec       X
295:       PetscInt  i,j,row
296:       PetscInt  xs,ys,xm
297:       PetscInt  ym
298:       PetscReal one,lambda,temp1,temp,hx,hy
299:       PetscScalar      xx(2)

301:       one    = 1.0
302:       lambda = 6.0
303:       hx     = one/(mx-1)
304:       hy     = one/(my-1)
305:       temp1  = lambda/(lambda + one)

307: !  Get a pointer to vector data.
308: !    - VecGetArray() returns a pointer to the data array.
309: !    - You MUST call VecRestoreArray() when you no longer need access to
310: !      the array.
311:        call VecGetArray(X,xx,idx,ierr)

313: !  Get local grid boundaries (for 2-dimensional DMDA):
314: !    xs, ys   - starting grid indices (no ghost points)
315: !    xm, ym   - widths of local grid (no ghost points)

317:        call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr)

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

321:       do 30 j=ys,ys+ym-1
322:         temp = (min(j,my-j-1))*hy
323:         do 40 i=xs,xs+xm-1
324:           row = i - xs + (j - ys)*xm + 1
325:           if (i .eq. 0 .or. j .eq. 0 .or. i .eq. mx-1 .or. j .eq. my-1) then
326:             xx(idx+row) = 0.0
327:             continue
328:           endif
329:           xx(idx+row) = temp1*sqrt(min((min(i,mx-i-1))*hx,temp))
330:  40     continue
331:  30   continue

333: !     Restore vector

335:        call VecRestoreArray(X,xx,idx,ierr)
336:        return
337:        end

339: ! -------------------------------------------------------------------
340: !
341: !   ComputeFunction - Evaluates nonlinear function, F(x).
342: !
343: !   Input Parameters:
344: !.  X - input vector
345: !
346: !   Output Parameter:
347: !.  F - function vector
348: !
349:       subroutine  ComputeFunction(X,F,ierr)
350:       use mymoduleex14f
351:       implicit none

353:       Vec              X,F
354:       PetscInt         gys,gxm,gym
355:       PetscOffset      idx,idf
356:       PetscErrorCode ierr
357:       PetscInt i,j,row,xs,ys,xm,ym,gxs
358:       PetscInt rowf
359:       PetscReal two,one,lambda,hx
360:       PetscReal hy,hxdhy,hydhx,sc
361:       PetscScalar      u,uxx,uyy,xx(2),ff(2)

363:       two    = 2.0
364:       one    = 1.0
365:       lambda = 6.0

367:       hx     = one/(mx-1)
368:       hy     = one/(my-1)
369:       sc     = hx*hy*lambda
370:       hxdhy  = hx/hy
371:       hydhx  = hy/hx

373: !  Scatter ghost points to local vector, using the 2-step process
374: !     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
375: !  By placing code between these two statements, computations can be
376: !  done while messages are in transition.
377: !
378:       call DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)
379:       call DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)

381: !  Get pointers to vector data

383:       call VecGetArray(localX,xx,idx,ierr)
384:       call VecGetArray(F,ff,idf,ierr)

386: !  Get local grid boundaries

388:       call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr)
389:       call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr)

391: !  Compute function over the locally owned part of the grid
392:       rowf = 0
393:       do 50 j=ys,ys+ym-1

395:         row  = (j - gys)*gxm + xs - gxs
396:         do 60 i=xs,xs+xm-1
397:           row  = row + 1
398:           rowf = rowf + 1

400:           if (i .eq. 0 .or. j .eq. 0 .or. i .eq. mx-1 .or. j .eq. my-1) then
401:             ff(idf+rowf) = xx(idx+row)
402:             goto 60
403:           endif
404:           u   = xx(idx+row)
405:           uxx = (two*u - xx(idx+row-1) - xx(idx+row+1))*hydhx
406:           uyy = (two*u - xx(idx+row-gxm) - xx(idx+row+gxm))*hxdhy
407:           ff(idf+rowf) = uxx + uyy - sc*exp(u)
408:  60     continue
409:  50   continue

411: !  Restore vectors

413:        call VecRestoreArray(localX,xx,idx,ierr)
414:        call VecRestoreArray(F,ff,idf,ierr)
415:        return
416:        end

418: ! -------------------------------------------------------------------
419: !
420: !   ComputeJacobian - Evaluates Jacobian matrix.
421: !
422: !   Input Parameters:
423: !   x - input vector
424: !
425: !   Output Parameters:
426: !   jac - Jacobian matrix
427: !   flag - flag indicating matrix structure
428: !
429: !   Notes:
430: !   Due to grid point reordering with DMDAs, we must always work
431: !   with the local grid points, and then transform them to the new
432: !   global numbering with the 'ltog' mapping
433: !   We cannot work directly with the global numbers for the original
434: !   uniprocessor grid!
435: !
436:       subroutine ComputeJacobian(X,jac,ierr)
437:       use mymoduleex14f
438:       implicit none

440:       Vec         X
441:       Mat         jac
442:       PetscInt     ltog(2)
443:       PetscOffset idltog,idx
444:       PetscErrorCode ierr
445:       PetscInt xs,ys,xm,ym
446:       PetscInt gxs,gys,gxm,gym
447:       PetscInt grow(1),i,j
448:       PetscInt row,ione
449:       PetscInt col(5),ifive
450:       PetscScalar two,one,lambda
451:       PetscScalar v(5),hx,hy,hxdhy
452:       PetscScalar hydhx,sc,xx(2)
453:       ISLocalToGlobalMapping ltogm

455:       ione   = 1
456:       ifive  = 5
457:       one    = 1.0
458:       two    = 2.0
459:       hx     = one/(mx-1)
460:       hy     = one/(my-1)
461:       sc     = hx*hy
462:       hxdhy  = hx/hy
463:       hydhx  = hy/hx
464:       lambda = 6.0

466: !  Scatter ghost points to local vector, using the 2-step process
467: !     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
468: !  By placing code between these two statements, computations can be
469: !  done while messages are in transition.

471:       call DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)
472:       call DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)

474: !  Get pointer to vector data

476:       call VecGetArray(localX,xx,idx,ierr)

478: !  Get local grid boundaries

480:       call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr)
481:       call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr)

483: !  Get the global node numbers for all local nodes, including ghost points

485:       call DMGetLocalToGlobalMapping(da,ltogm,ierr)
486:       call ISLocalToGlobalMappingGetIndices(ltogm,ltog,idltog,ierr)

488: !  Compute entries for the locally owned part of the Jacobian.
489: !   - Currently, all PETSc parallel matrix formats are partitioned by
490: !     contiguous chunks of rows across the processors. The 'grow'
491: !     parameter computed below specifies the global row number
492: !     corresponding to each local grid point.
493: !   - Each processor needs to insert only elements that it owns
494: !     locally (but any non-local elements will be sent to the
495: !     appropriate processor during matrix assembly).
496: !   - Always specify global row and columns of matrix entries.
497: !   - Here, we set all entries for a particular row at once.

499:       do 10 j=ys,ys+ym-1
500:         row = (j - gys)*gxm + xs - gxs
501:         do 20 i=xs,xs+xm-1
502:           row = row + 1
503:           grow(1) = ltog(idltog+row)
504:           if (i .eq. 0 .or. j .eq. 0 .or. i .eq. (mx-1) .or. j .eq. (my-1)) then
505:              call MatSetValues(jac,ione,grow,ione,grow,one,INSERT_VALUES,ierr)
506:              go to 20
507:           endif
508:           v(1)   = -hxdhy
509:           col(1) = ltog(idltog+row - gxm)
510:           v(2)   = -hydhx
511:           col(2) = ltog(idltog+row - 1)
512:           v(3)   = two*(hydhx + hxdhy) - sc*lambda*exp(xx(idx+row))
513:           col(3) = grow(1)
514:           v(4)   = -hydhx
515:           col(4) = ltog(idltog+row + 1)
516:           v(5)   = -hxdhy
517:           col(5) = ltog(idltog+row + gxm)
518:           call MatSetValues(jac,ione,grow,ifive,col,v,INSERT_VALUES,ierr)
519:  20     continue
520:  10   continue

522:       call ISLocalToGlobalMappingRestoreIndices(ltogm,ltog,idltog,ierr)

524: !  Assemble matrix, using the 2-step process:
525: !    MatAssemblyBegin(), MatAssemblyEnd().
526: !  By placing code between these two statements, computations can be
527: !  done while messages are in transition.

529:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
530:       call VecRestoreArray(localX,xx,idx,ierr)
531:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
532:       return
533:       end


536: ! -------------------------------------------------------------------
537: !
538: !   MyMult - user provided matrix multiply
539: !
540: !   Input Parameters:
541: !.  X - input vector
542: !
543: !   Output Parameter:
544: !.  F - function vector
545: !
546:       subroutine  MyMult(J,X,F,ierr)
547:       use mymoduleex14f
548:       implicit none
549: 
550:       Mat     J
551:       Vec     X,F
552:       PetscErrorCode ierr
553: !
554: !       Here we use the actual formed matrix B; users would
555: !     instead write their own matrix vector product routine
556: !
557:       call MatMult(B,X,F,ierr)
558:       return
559:       end

561: !/*TEST
562: !
563: !   test:
564: !      args: -no_output -ksp_gmres_cgs_refinement_type refine_always
565: !      output_file: output/ex14_1.out
566: !      requires: !single
567: !
568: !TEST*/