Actual source code: ex15f.F

petsc-master 2016-09-28
Report Typos and Errors
  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: !
  7: !/*T
  8: !   Concepts: KSP^basic parallel example
  9: !   Concepts: PC^setting a user-defined shell preconditioner
 10: !   Processors: n
 11: !T*/
 12: !
 13: !  -------------------------------------------------------------------------

 15:       program main
 16:       implicit none

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

 26:  #include <petsc/finclude/petscsys.h>
 27:  #include <petsc/finclude/petscvec.h>
 28:  #include <petsc/finclude/petscmat.h>
 29:  #include <petsc/finclude/petscpc.h>
 30:  #include <petsc/finclude/petscksp.h>

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

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

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

 59:       external SampleShellPCSetUp, SampleShellPCApply
 60:       external  SampleShellPCDestroy

 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:       if (ierr .ne. 0) then
 72:         print*,'Unable to initialize PETSc'
 73:         stop
 74:       endif
 75:       one     = 1.0
 76:       neg_one = -1.0
 77:       i1 = 1
 78:       m       = 8
 79:       n       = 7
 80:       five    = 5
 81:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,     &
 82:      &                        '-m',m,flg,ierr)
 83:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,     &
 84:      &                        '-n',n,flg,ierr)
 85:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 87: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88: !      Compute the matrix and right-hand-side vector that define
 89: !      the linear system, Ax = b.
 90: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 92: !  Create parallel matrix, specifying only its global dimensions.
 93: !  When using MatCreate(), the matrix format can be specified at
 94: !  runtime. Also, the parallel partitioning of the matrix is
 95: !  determined by PETSc at runtime.

 97:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 98:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,ierr)
 99:       call MatSetType(A, MATAIJ,ierr)
100:       call MatSetFromOptions(A,ierr)
101:       call MatMPIAIJSetPreallocation(A,five,PETSC_NULL_INTEGER,five,            &
102:      &                     PETSC_NULL_INTEGER,ierr)
103:       call MatSeqAIJSetPreallocation(A,five,PETSC_NULL_INTEGER,ierr)

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

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

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

119:       do 10, II=Istart,Iend-1
120:         v = -1.0
121:         i = II/n
122:         j = II - i*n
123:         if (i.gt.0) then
124:           JJ = II - n
125:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
126:         endif
127:         if (i.lt.m-1) then
128:           JJ = II + n
129:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
130:         endif
131:         if (j.gt.0) then
132:           JJ = II - 1
133:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
134:         endif
135:         if (j.lt.n-1) then
136:           JJ = II + 1
137:           call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
138:         endif
139:         v = 4.0
140:         call  MatSetValues(A,i1,II,i1,II,v,ADD_VALUES,ierr)
141:  10   continue

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

148:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
149:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

151: !  Create parallel vectors.
152: !   - Here, the parallel partitioning of the vector is determined by
153: !     PETSc at runtime.  We could also specify the local dimensions
154: !     if desired -- or use the more general routine VecCreate().
155: !   - When solving a linear system, the vectors and matrices MUST
156: !     be partitioned accordingly.  PETSc automatically generates
157: !     appropriately partitioned matrices and vectors when MatCreate()
158: !     and VecCreate() are used with the same communicator.
159: !   - Note: We form 1 vector from scratch and then duplicate as needed.

161:       call VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,m*n,u,ierr)
162:       call VecDuplicate(u,b,ierr)
163:       call VecDuplicate(b,x,ierr)

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

167:       call VecSet(u,one,ierr)
168:       call MatMult(A,u,b,ierr)

170: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171: !         Create the linear solver and set various options
172: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

174: !  Create linear solver context

176:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

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

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

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

188:       call KSPGetPC(ksp,pc,ierr)
189:       tol = 1.e-7
190:       call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_REAL,                       &
191:      &     PETSC_DEFAULT_REAL,PETSC_DEFAULT_INTEGER,ierr)

193: !
194: !  Set a user-defined shell preconditioner if desired
195: !
196:       call PetscOptionsHasName(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,        &
197:      &                         '-user_defined_pc',user_defined_pc,ierr)

199:       if (user_defined_pc) then

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

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

207: !  (Optional) Do any setup required for the preconditioner
208:          call PCShellSetSetUp(pc,SampleShellPCSetUp,ierr)

210: !  (Optional) Frees any objects we created for the preconditioner
211:          call PCShellSetDestroy(pc,SampleShellPCDestroy,ierr)

213:       else
214:          call PCSetType(pc,PCJACOBI,ierr)
215:       endif

217: !  Set runtime options, e.g.,
218: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
219: !  These options will override those specified above as long as
220: !  KSPSetFromOptions() is called _after_ any other customization
221: !  routines.

223:       call KSPSetFromOptions(ksp,ierr)

225: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226: !                      Solve the linear system
227: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

231: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232: !                     Check solution and clean up
233: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

235: !  Check the error

237:       call VecAXPY(x,neg_one,u,ierr)
238:       call VecNorm(x,NORM_2,norm,ierr)
239:       call KSPGetIterationNumber(ksp,its,ierr)

241:       if (rank .eq. 0) then
242:         if (norm .gt. 1.e-12) then
243:            write(6,100) norm,its
244:         else
245:            write(6,110) its
246:         endif
247:       endif
248:   100 format('Norm of error ',1pe11.4,' iterations ',i5)
249:   110 format('Norm of error < 1.e-12,iterations ',i5)

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

254:       call KSPDestroy(ksp,ierr)
255:       call VecDestroy(u,ierr)
256:       call VecDestroy(x,ierr)
257:       call VecDestroy(b,ierr)
258:       call MatDestroy(A,ierr)

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

262:       call PetscFinalize(ierr)
263:       end

265: !/***********************************************************************/
266: !/*          Routines for a user-defined shell preconditioner           */
267: !/***********************************************************************/

269: !
270: !   SampleShellPCSetUp - This routine sets up a user-defined
271: !   preconditioner context.
272: !
273: !   Input Parameters:
274: !   pc - preconditioner object
275: !
276: !   Output Parameter:
277: !   ierr  - error code (nonzero if error has been detected)
278: !
279: !   Notes:
280: !   In this example, we define the shell preconditioner to be Jacobi
281: !   method.  Thus, here we create a work vector for storing the reciprocal
282: !   of the diagonal of the preconditioner matrix; this vector is then
283: !   used within the routine SampleShellPCApply().
284: !
285:       subroutine SampleShellPCSetUp(pc,ierr)

287:       implicit none

289:  #include <petsc/finclude/petscsys.h>
290:  #include <petsc/finclude/petscvec.h>
291:  #include <petsc/finclude/petscmat.h>
292:       PC      pc

294:       Mat     pmat
295:       integer ierr

297: !  Common block to store data for user-provided preconditioner
298: !  Normally we would recommend storing all the work data (like diag) in
299: !  the context set with PCShellSetContext()

301:       common /myshellpc/ diag
302:       Vec    diag

304:       call PCGetOperators(pc,PETSC_NULL_OBJECT,pmat,ierr)
305:       call MatCreateVecs(pmat,diag,PETSC_NULL_OBJECT,ierr)
306:       call MatGetDiagonal(pmat,diag,ierr)
307:       call VecReciprocal(diag,ierr)

309:       end

311: ! -------------------------------------------------------------------
312: !
313: !   SampleShellPCApply - This routine demonstrates the use of a
314: !   user-provided preconditioner.
315: !
316: !   Input Parameters:
317: !   pc - preconditioner object
318: !   x - input vector
319: !
320: !   Output Parameters:
321: !   y - preconditioned vector
322: !   ierr  - error code (nonzero if error has been detected)
323: !
324: !   Notes:
325: !   This code implements the Jacobi preconditioner, merely as an
326: !   example of working with a PCSHELL.  Note that the Jacobi method
327: !   is already provided within PETSc.
328: !
329:       subroutine SampleShellPCApply(pc,x,y,ierr)

331:       implicit none

333:  #include <petsc/finclude/petscsys.h>
334:  #include <petsc/finclude/petscvec.h>

336:       PC      pc
337:       Vec     x,y
338:       integer ierr

340: !  Common block to store data for user-provided preconditioner
341:       common /myshellpc/ diag
342:       Vec    diag

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

346:       end

348: !/***********************************************************************/
349: !/*          Routines for a user-defined shell preconditioner           */
350: !/***********************************************************************/

352: !
353: !   SampleShellPCDestroy - This routine destroys (frees the memory of) any
354: !      objects we made for the preconditioner
355: !
356: !   Input Parameters:
357: !   pc - for this example we use the actual PC as our shell context
358: !
359: !   Output Parameter:
360: !   ierr  - error code (nonzero if error has been detected)
361: !

363:       subroutine SampleShellPCDestroy(pc,ierr)

365:       implicit none

367:  #include <petsc/finclude/petscsys.h>
368:  #include <petsc/finclude/petscvec.h>
369:  #include <petsc/finclude/petscmat.h>
370:       PC      pc
371:       PetscErrorCode ierr

373: !  Common block to store data for user-provided preconditioner
374: !  Normally we would recommend storing all the work data (like diag) in
375: !  the context set with PCShellSetContext()

377:       common /myshellpc/ diag
378:       Vec    diag

380:       call VecDestroy(diag,ierr)

382:       end