Actual source code: ex11f.F

petsc-3.4.5 2014-06-29
  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 <finclude/petscsys.h>
 43: #include <finclude/petscvec.h>
 44: #include <finclude/petscmat.h>
 45: #include <finclude/petscpc.h>
 46: #include <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:       double precision 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 !defined(PETSC_USE_COMPLEX)
 81:       write(6,*) "This example requires complex numbers."
 82:       goto 200
 83: #endif

 85:       none   = -1.0
 86:       n      = 6
 87:       sigma1 = 100.0
 88:       czero  = 0.0
 89:       cone   = PETSC_i
 90:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 91:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-sigma1',sigma1,      &
 92:      &                       flg,ierr)
 93:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
 94:       dim    = n*n

 96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97: !      Compute the matrix and right-hand-side vector that define
 98: !      the linear system, Ax = b.
 99: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

101: !  Create parallel matrix, specifying only its global dimensions.
102: !  When using MatCreate(), the matrix format can be specified at
103: !  runtime. Also, the parallel partitioning of the matrix is
104: !  determined by PETSc at runtime.

106:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
107:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim,ierr)
108:       call MatSetFromOptions(A,ierr)
109:       call MatSetUp(A,ierr)

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

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

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

123:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-norandom',        &
124:      &     flg,ierr)
125:       if (flg) then
126:          use_random = .false.
127:          sigma2 = 10.0*PETSC_i
128:       else
129:          use_random = .true.
130:          call PetscRandomCreate(PETSC_COMM_WORLD,                       &
131:      &        rctx,ierr)
132:          call PetscRandomSetFromOptions(rctx,ierr)
133:          call PetscRandomSetInterval(rctx,czero,cone,ierr)
134:       endif
135:       h2 = 1.0/((n+1)*(n+1))

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

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

170:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
171:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

173: !  Create parallel vectors.
174: !   - Here, the parallel partitioning of the vector is determined by
175: !     PETSc at runtime.  We could also specify the local dimensions
176: !     if desired.
177: !   - Note: We form 1 vector from scratch and then duplicate as needed.

179:       call VecCreate(PETSC_COMM_WORLD,u,ierr)
180:       call VecSetSizes(u,PETSC_DECIDE,dim,ierr)
181:       call VecSetFromOptions(u,ierr)
182:       call VecDuplicate(u,b,ierr)
183:       call VecDuplicate(b,x,ierr)

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

187:       if (use_random) then
188:          call PetscRandomCreate(PETSC_COMM_WORLD,rctx,ierr)
189:          call PetscRandomSetFromOptions(rctx,ierr)
190:          call VecSetRandom(u,rctx,ierr)
191:       else
192:          pfive = 0.5
193:          call VecSet(u,pfive,ierr)
194:       endif
195:       call MatMult(A,u,b,ierr)

197: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: !         Create the linear solver and set various options
199: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

201: !  Create linear solver context

203:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

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

208:       call KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN,ierr)

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

213:       call KSPSetFromOptions(ksp,ierr)

215: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216: !                      Solve the linear system
217: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

221: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222: !                     Check solution and clean up
223: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

225: !  Check the error

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

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

243:       if (use_random) call PetscRandomDestroy(rctx,ierr)
244:       call KSPDestroy(ksp,ierr)
245:       call VecDestroy(u,ierr)
246:       call VecDestroy(x,ierr)
247:       call VecDestroy(b,ierr)
248:       call MatDestroy(A,ierr)

250: #if !defined(PETSC_USE_COMPLEX)
251:  200  continue
252: #endif
253:       call PetscFinalize(ierr)
254:       end