Actual source code: ex15f.F

  1: !
  2: !   Solves a linear system in parallel with KSP.  Also indicates
  3: !   use of a user-provided preconditioner.  Input parameters include:
  4: !      -user_defined_pc : Activate a user-defined preconditioner
  5: !
  6: !  Program usage: mpirun ex15f [-help] [all PETSc options]
  7: !
  8: !/*T
  9: !   Concepts: KSP^basic parallel example
 10: !   Concepts: PC^setting a user-defined shell preconditioner
 11: !   Processors: n
 12: !T*/
 13: !
 14: !  -------------------------------------------------------------------------

 16:       program main
 17:       implicit none

 19: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: !                    Include files
 21: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: !
 23: !     petsc.h  - base PETSc routines      petscvec.h - vectors
 24: !     petscsys.h    - system routines          petscmat.h - matrices
 25: !     petscksp.h    - Krylov subspace methods  petscpc.h  - preconditioners

 27:  #include include/finclude/petsc.h
 28:  #include include/finclude/petscvec.h
 29:  #include include/finclude/petscmat.h
 30:  #include include/finclude/petscpc.h
 31:  #include include/finclude/petscksp.h

 33: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 34: !                   Variable declarations
 35: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 36: !
 37: !  Variables:
 38: !     ksp     - linear solver context
 39: !     ksp      - Krylov subspace method context
 40: !     pc       - preconditioner context
 41: !     x, b, u  - approx solution, right-hand-side, exact solution vectors
 42: !     A        - matrix that defines linear system
 43: !     its      - iterations for convergence
 44: !     norm     - norm of solution error

 46:       Vec              x,b,u
 47:       Mat              A
 48:       PC               pc
 49:       KSP              ksp
 50:       PetscScalar      v,one,neg_one
 51:       double precision norm,tol
 52:       PetscErrorCode ierr
 53:       PetscInt   i,j,II,JJ,Istart
 54:       PetscInt   Iend,m,n,i1,its
 55:       PetscMPIInt rank
 56:       PetscTruth user_defined_pc,flg

 58: !  Note: Any user-defined Fortran routines MUST be declared as external.

 60:       external SampleShellPCSetUp, SampleShellPCApply

 62: !  Common block to store data for user-provided preconditioner
 63:       common /myshellpc/ diag
 64:       Vec    diag

 66: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67: !                 Beginning of program
 68: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 70:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 71:       one     = 1.0
 72:       neg_one = -1.0
 73:       i1 = 1
 74:       m       = 8
 75:       n       = 7
 76:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-m',m,flg,ierr)
 77:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
 78:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 80: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 81: !      Compute the matrix and right-hand-side vector that define
 82: !      the linear system, Ax = b.
 83: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 85: !  Create parallel matrix, specifying only its global dimensions.
 86: !  When using MatCreate(), the matrix format can be specified at
 87: !  runtime. Also, the parallel partitioning of the matrix is
 88: !  determined by PETSc at runtime.

 90:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 91:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,ierr)
 92:       call MatSetFromOptions(A,ierr)

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

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

100: !  Set matrix elements for the 2-D, five-point stencil in parallel.
101: !   - Each processor needs to insert only elements that it owns
102: !     locally (but any non-local elements will be sent to the
103: !     appropriate processor during matrix assembly).
104: !   - Always specify global row and columns of matrix entries.
105: !   - Note that MatSetValues() uses 0-based row and column numbers
106: !     in Fortran as well as in C.

108:       do 10, II=Istart,Iend-1
109:         v = -1.0
110:         i = II/n
111:         j = II - i*n
112:         if (i.gt.0) then
113:           JJ = II - n
114:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
115:         endif
116:         if (i.lt.m-1) then
117:           JJ = II + n
118:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
119:         endif
120:         if (j.gt.0) then
121:           JJ = II - 1
122:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
123:         endif
124:         if (j.lt.n-1) then
125:           JJ = II + 1
126:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
127:         endif
128:         v = 4.0
129:         call  MatSetValues(A,i1,II,i1,II,v,ADD_VALUES,ierr)
130:  10   continue

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

137:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
138:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

140: !  Create parallel vectors.
141: !   - Here, the parallel partitioning of the vector is determined by
142: !     PETSc at runtime.  We could also specify the local dimensions
143: !     if desired -- or use the more general routine VecCreate().
144: !   - When solving a linear system, the vectors and matrices MUST
145: !     be partitioned accordingly.  PETSc automatically generates
146: !     appropriately partitioned matrices and vectors when MatCreate()
147: !     and VecCreate() are used with the same communicator.
148: !   - Note: We form 1 vector from scratch and then duplicate as needed.

150:       call VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,m*n,u,ierr)
151:       call VecDuplicate(u,b,ierr)
152:       call VecDuplicate(b,x,ierr)

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

156:       call VecSet(u,one,ierr)
157:       call MatMult(A,u,b,ierr)

159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: !         Create the linear solver and set various options
161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

163: !  Create linear solver context

165:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

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

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

172: !  Set linear solver defaults for this problem (optional).
173: !   - By extracting the KSP and PC contexts from the KSP context,
174: !     we can then directly directly call any KSP and PC routines
175: !     to set various options.

177:       call KSPGetPC(ksp,pc,ierr)
178:       tol = 1.e-7
179:       call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_DOUBLE_PRECISION,     &
180:      &     PETSC_DEFAULT_DOUBLE_PRECISION,PETSC_DEFAULT_INTEGER,ierr)

182: !
183: !  Set a user-defined shell preconditioner if desired
184: !
185:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-user_defined_pc',      &
186:      &                    user_defined_pc,ierr)

188:       if (user_defined_pc .eq. 1) then

190: !  (Required) Indicate to PETSc that we are using a shell preconditioner
191:          call PCSetType(pc,PCSHELL,ierr)

193: !  (Required) Set the user-defined routine for applying the preconditioner
194:          call PCShellSetApply(pc,SampleShellPCApply,ierr)

196: !  (Optional) Do any setup required for the preconditioner
197:          call SampleShellPCSetUp(A,x,ierr)

199:       else
200:          call PCSetType(pc,PCJACOBI,ierr)
201:       endif

203: !  Set runtime options, e.g.,
204: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
205: !  These options will override those specified above as long as
206: !  KSPSetFromOptions() is called _after_ any other customization
207: !  routines.

209:       call KSPSetFromOptions(ksp,ierr)

211: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212: !                      Solve the linear system
213: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

217: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218: !                     Check solution and clean up
219: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

221: !  Check the error

223:       call VecAXPY(x,neg_one,u,ierr)
224:       call VecNorm(x,NORM_2,norm,ierr)
225:       call KSPGetIterationNumber(ksp,its,ierr)

227:       if (rank .eq. 0) then
228:         if (norm .gt. 1.e-12) then
229:            write(6,100) norm,its
230:         else
231:            write(6,110) its
232:         endif
233:       endif
234:   100 format('Norm of error ',1pe10.4,' iterations ',i5)
235:   110 format('Norm of error < 1.e-12,iterations ',i5)

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

240:       call KSPDestroy(ksp,ierr)
241:       call VecDestroy(u,ierr)
242:       call VecDestroy(x,ierr)
243:       call VecDestroy(b,ierr)
244:       call MatDestroy(A,ierr)
245:       if (user_defined_pc .eq. 1) then
246:          call VecDestroy(diag,ierr)
247:       endif


250: !  Always call PetscFinalize() before exiting a program.

252:       call PetscFinalize(ierr)
253:       end

255: !/***********************************************************************/
256: !/*          Routines for a user-defined shell preconditioner           */
257: !/***********************************************************************/

259: !
260: !   SampleShellPCSetUp - This routine sets up a user-defined
261: !   preconditioner context.
262: !
263: !   Input Parameters:
264: !   pmat  - preconditioner matrix
265: !   x     - vector
266: !
267: !   Output Parameter:
268: !   ierr  - error code (nonzero if error has been detected)
269: !
270: !   Notes:
271: !   In this example, we define the shell preconditioner to be Jacobi
272: !   method.  Thus, here we create a work vector for storing the reciprocal
273: !   of the diagonal of the preconditioner matrix; this vector is then
274: !   used within the routine SampleShellPCApply().
275: !
276:       subroutine SampleShellPCSetUp(pmat,x,ierr)

278:       implicit none

280:  #include include/finclude/petsc.h
281:  #include include/finclude/petscvec.h
282:  #include include/finclude/petscmat.h

284:       Vec     x
285:       Mat     pmat
286:       integer ierr

288: !  Common block to store data for user-provided preconditioner
289:       common /myshellpc/ diag
290:       Vec    diag

292:       call VecDuplicate(x,diag,ierr)
293:       call MatGetDiagonal(pmat,diag,ierr)
294:       call VecReciprocal(diag,ierr)

296:       end

298: ! -------------------------------------------------------------------
299: !
300: !   SampleShellPCApply - This routine demonstrates the use of a
301: !   user-provided preconditioner.
302: !
303: !   Input Parameters:
304: !   dummy - optional user-defined context, not used here
305: !   x - input vector
306: !
307: !   Output Parameters:
308: !   y - preconditioned vector
309: !   ierr  - error code (nonzero if error has been detected)
310: !
311: !   Notes:
312: !   This code implements the Jacobi preconditioner, merely as an
313: !   example of working with a PCSHELL.  Note that the Jacobi method
314: !   is already provided within PETSc.
315: !
316:       subroutine SampleShellPCApply(dummy,x,y,ierr)

318:       implicit none

320:  #include include/finclude/petsc.h
321:  #include include/finclude/petscvec.h

323:       Vec     x,y
324:       integer dummy,ierr

326: !  Common block to store data for user-provided preconditioner
327:       common /myshellpc/ diag
328:       Vec    diag

330:       call VecPointwiseMult(y,x,diag,ierr)

332:       end