Actual source code: ex11f.F

petsc-master 2016-12-04
Report Typos and Errors
  1: !
  2: !  Description: Solves a complex linear system in parallel with KSP (Fortran code).
  3: !
  4: !/*T
  5: !  Concepts: KSP^solving a Helmholtz equation
  6: !  Concepts: complex numbers
  7: !  Processors: n
  8: !T*/
  9: !
 10: !  The model problem:
 11: !     Solve Helmholtz equation on the unit square: (0,1) x (0,1)
 12: !          -delta u - sigma1*u + i*sigma2*u = f,
 13: !           where delta = Laplace operator
 14: !     Dirichlet b.c.'s on all sides
 15: !     Use the 2-D, five-point finite difference stencil.
 16: !
 17: !     Compiling the code:
 18: !      This code uses the complex numbers version of PETSc, so configure
 19: !      must be run to enable this
 20: !
 21: !
 22: ! -----------------------------------------------------------------------

 24:       program main
 25:       implicit none

 27: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 28: !                    Include files
 29: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 30: !
 31: !  The following include statements are required for KSP Fortran programs:
 32: !     petscsys.h       - base PETSc routines
 33: !     petscvec.h    - vectors
 34: !     petscmat.h    - matrices
 35: !     petscpc.h     - preconditioners
 36: !     petscksp.h    - Krylov subspace methods
 37: !  Additional include statements may be needed if using other PETSc
 38: !  routines in a Fortran program, e.g.,
 39: !     petscviewer.h - viewers
 40: !     petscis.h     - index sets
 41: !
 42:  #include <petsc/finclude/petscsys.h>
 43:  #include <petsc/finclude/petscvec.h>
 44:  #include <petsc/finclude/petscmat.h>
 45:  #include <petsc/finclude/petscpc.h>
 46:  #include <petsc/finclude/petscksp.h>
 47: !
 48: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49: !                   Variable declarations
 50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51: !
 52: !  Variables:
 53: !     ksp     - linear solver context
 54: !     x, b, u  - approx solution, right-hand-side, exact solution vectors
 55: !     A        - matrix that defines linear system
 56: !     its      - iterations for convergence
 57: !     norm     - norm of error in solution
 58: !     rctx     - random number context
 59: !

 61:       KSP             ksp
 62:       Mat              A
 63:       Vec              x,b,u
 64:       PetscRandom      rctx
 65:       PetscReal norm,h2,sigma1
 66:       PetscScalar  none,sigma2,v,pfive,czero
 67:       PetscScalar  cone
 68:       PetscInt dim,its,n,Istart
 69:       PetscInt Iend,i,j,II,JJ,one
 70:       PetscErrorCode ierr
 71:       PetscMPIInt rank
 72:       PetscBool  flg
 73:       logical          use_random

 75: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76: !                 Beginning of program
 77: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 79:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 80:       if (ierr .ne. 0) then
 81:         print*,'Unable to initialize PETSc'
 82:         stop
 83:       endif
 84: #if !defined(PETSC_USE_COMPLEX)
 85:       write(6,*) 'This example requires complex numbers.'
 86:       goto 200
 87: #endif

 89:       none   = -1.0
 90:       n      = 6
 91:       sigma1 = 100.0
 92:       czero  = 0.0
 93:       cone   = PETSC_i
 94:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 95:       call PetscOptionsGetReal(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,     &
 96:      &                         '-sigma1',sigma1,flg,ierr)
 97:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,      &
 98:      &                        '-n',n,flg,ierr)
 99:       dim    = n*n

101: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: !      Compute the matrix and right-hand-side vector that define
103: !      the linear system, Ax = b.
104: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

106: !  Create parallel matrix, specifying only its global dimensions.
107: !  When using MatCreate(), the matrix format can be specified at
108: !  runtime. Also, the parallel partitioning of the matrix is
109: !  determined by PETSc at runtime.

111:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
112:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim,ierr)
113:       call MatSetFromOptions(A,ierr)
114:       call MatSetUp(A,ierr)

116: !  Currently, all PETSc parallel matrix formats are partitioned by
117: !  contiguous chunks of rows across the processors.  Determine which
118: !  rows of the matrix are locally owned.

120:       call MatGetOwnershipRange(A,Istart,Iend,ierr)

122: !  Set matrix elements in parallel.
123: !   - Each processor needs to insert only elements that it owns
124: !     locally (but any non-local elements will be sent to the
125: !     appropriate processor during matrix assembly).
126: !   - Always specify global rows and columns of matrix entries.

128:       call PetscOptionsHasName(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,   &
129:      &                         '-norandom',flg,ierr)
130:       if (flg) then
131:          use_random = .false.
132:          sigma2 = 10.0*PETSC_i
133:       else
134:          use_random = .true.
135:          call PetscRandomCreate(PETSC_COMM_WORLD,                       &
136:      &        rctx,ierr)
137:          call PetscRandomSetFromOptions(rctx,ierr)
138:          call PetscRandomSetInterval(rctx,czero,cone,ierr)
139:       endif
140:       h2 = 1.0/real((n+1)*(n+1))

142:       one = 1
143:       do 10, II=Istart,Iend-1
144:         v = -1.0
145:         i = II/n
146:         j = II - i*n
147:         if (i.gt.0) then
148:           JJ = II - n
149:           call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
150:         endif
151:         if (i.lt.n-1) then
152:           JJ = II + n
153:           call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
154:         endif
155:         if (j.gt.0) then
156:           JJ = II - 1
157:           call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
158:         endif
159:         if (j.lt.n-1) then
160:           JJ = II + 1
161:           call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
162:         endif
163:         if (use_random) call PetscRandomGetValue(rctx,                          &
164:      &                        sigma2,ierr)
165:         v = 4.0 - sigma1*h2 + sigma2*h2
166:         call  MatSetValues(A,one,II,one,II,v,ADD_VALUES,ierr)
167:  10   continue
168:       if (use_random) call PetscRandomDestroy(rctx,ierr)

170: !  Assemble matrix, using the 2-step process:
171: !       MatAssemblyBegin(), MatAssemblyEnd()
172: !  Computations can be done while messages are in transition
173: !  by placing code between these two statements.

175:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
176:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

178: !  Create parallel vectors.
179: !   - Here, the parallel partitioning of the vector is determined by
180: !     PETSc at runtime.  We could also specify the local dimensions
181: !     if desired.
182: !   - Note: We form 1 vector from scratch and then duplicate as needed.

184:       call VecCreate(PETSC_COMM_WORLD,u,ierr)
185:       call VecSetSizes(u,PETSC_DECIDE,dim,ierr)
186:       call VecSetFromOptions(u,ierr)
187:       call VecDuplicate(u,b,ierr)
188:       call VecDuplicate(b,x,ierr)

190: !  Set exact solution; then compute right-hand-side vector.

192:       if (use_random) then
193:          call PetscRandomCreate(PETSC_COMM_WORLD,rctx,ierr)
194:          call PetscRandomSetFromOptions(rctx,ierr)
195:          call VecSetRandom(u,rctx,ierr)
196:       else
197:          pfive = 0.5
198:          call VecSet(u,pfive,ierr)
199:       endif
200:       call MatMult(A,u,b,ierr)

202: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: !         Create the linear solver and set various options
204: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

206: !  Create linear solver context

208:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

210: !  Set operators. Here the matrix that defines the linear system
211: !  also serves as the preconditioning matrix.

213:       call KSPSetOperators(ksp,A,A,ierr)

215: !  Set runtime options, e.g.,
216: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>

218:       call KSPSetFromOptions(ksp,ierr)

220: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221: !                      Solve the linear system
222: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

224:       call KSPSolve(ksp,b,x,ierr)

226: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227: !                     Check solution and clean up
228: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

230: !  Check the error

232:       call VecAXPY(x,none,u,ierr)
233:       call VecNorm(x,NORM_2,norm,ierr)
234:       call KSPGetIterationNumber(ksp,its,ierr)
235:       if (rank .eq. 0) then
236:         if (norm .gt. 1.e-12) then
237:            write(6,100) norm,its
238:         else
239:            write(6,110) its
240:         endif
241:       endif
242:   100 format('Norm of error ',e11.4,',iterations ',i5)
243:   110 format('Norm of error < 1.e-12,iterations ',i5)

245: !  Free work space.  All PETSc objects should be destroyed when they
246: !  are no longer needed.

248:       if (use_random) call PetscRandomDestroy(rctx,ierr)
249:       call KSPDestroy(ksp,ierr)
250:       call VecDestroy(u,ierr)
251:       call VecDestroy(x,ierr)
252:       call VecDestroy(b,ierr)
253:       call MatDestroy(A,ierr)

255: #if !defined(PETSC_USE_COMPLEX)
256:  200  continue
257: #endif
258:       call PetscFinalize(ierr)
259:       end