Actual source code: ex15f.F90

  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: !
  8: !     Module contains diag needed by shell preconditioner
  9: !
 10: #include <petsc/finclude/petscksp.h>
 11: module ex15fmodule
 12:   use petscksp
 13:   implicit none
 14:   Vec diag

 16: contains

 18: !/***********************************************************************/
 19: !/*          Routines for a user-defined shell preconditioner           */
 20: !/***********************************************************************/

 22: !
 23: !   SampleShellPCSetUp - This routine sets up a user-defined
 24: !   preconditioner context.
 25: !
 26: !   Input Parameters:
 27: !   pc - preconditioner object
 28: !
 29: !   Output Parameter:
 30: !   ierr  - error code (nonzero if error has been detected)
 31: !
 32: !   Notes:
 33: !   In this example, we define the shell preconditioner to be Jacobi
 34: !   method.  Thus, here we create a work vector for storing the reciprocal
 35: !   of the diagonal of the matrix; this vector is then
 36: !   used within the routine SampleShellPCApply().
 37: !
 38:   subroutine SampleShellPCSetUp(pc, ierr)

 40:     PC pc
 41:     Mat pmat
 42:     PetscErrorCode, intent(out) :: ierr

 44:     PetscCallA(PCGetOperators(pc, PETSC_NULL_MAT, pmat, ierr))
 45:     PetscCallA(MatCreateVecs(pmat, diag, PETSC_NULL_VEC, ierr))
 46:     PetscCallA(MatGetDiagonal(pmat, diag, ierr))
 47:     PetscCallA(VecReciprocal(diag, ierr))

 49:   end

 51: ! -------------------------------------------------------------------
 52: !
 53: !   SampleShellPCApply - This routine demonstrates the use of a
 54: !   user-provided preconditioner.
 55: !
 56: !   Input Parameters:
 57: !   pc - preconditioner object
 58: !   x - input vector
 59: !
 60: !   Output Parameters:
 61: !   y - preconditioned vector
 62: !   ierr  - error code (nonzero if error has been detected)
 63: !
 64: !   Notes:
 65: !   This code implements the Jacobi preconditioner, merely as an
 66: !   example of working with a PCSHELL.  Note that the Jacobi method
 67: !   is already provided within PETSc.
 68: !
 69:   subroutine SampleShellPCApply(pc, x, y, ierr)
 70:     PC pc
 71:     Vec x, y
 72:     PetscErrorCode, intent(out) :: ierr

 74:     PetscCallA(VecPointwiseMult(y, x, diag, ierr))

 76:   end

 78: !/***********************************************************************/
 79: !/*          Routines for a user-defined shell preconditioner           */
 80: !/***********************************************************************/

 82: !
 83: !   SampleShellPCDestroy - This routine destroys (frees the memory of) any
 84: !      objects we made for the preconditioner
 85: !
 86: !   Input Parameters:
 87: !   pc - for this example we use the actual PC as our shell context
 88: !
 89: !   Output Parameter:
 90: !   ierr  - error code (nonzero if error has been detected)
 91: !

 93:   subroutine SampleShellPCDestroy(pc, ierr)
 94:     PC pc
 95:     PetscErrorCode, intent(out) :: ierr

 97: !  Normally we would recommend storing all the work data (like diag) in
 98: !  the context set with PCShellSetContext()

100:     PetscCallA(VecDestroy(diag, ierr))

102:   end

104: end module

106: program main
107:   use ex15fmodule
108:   implicit none

110: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: !                   Variable declarations
112: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: !
114: !  Variables:
115: !     ksp     - linear solver context
116: !     ksp      - Krylov subspace method context
117: !     pc       - preconditioner context
118: !     x, b, u  - approx solution, right-hand side, exact solution vectors
119: !     A        - matrix that defines linear system
120: !     its      - iterations for convergence
121: !     norm     - norm of solution error

123:   Vec x, b, u
124:   Mat A
125:   PC pc
126:   KSP ksp
127:   PetscScalar v
128:   PetscScalar, parameter :: one = 1.0, neg_one = -1.0
129:   PetscReal norm
130:   PetscReal, parameter :: tol = 1.e-7
131:   PetscErrorCode ierr
132:   PetscInt i, j, II, JJ, Istart, Iend, its
133:   PetscInt m, n
134:   PetscMPIInt rank
135:   PetscBool user_defined_pc, flg

137: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: !                 Beginning of program
139: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

141:   PetscCallA(PetscInitialize(ierr))

143:   m = 8
144:   PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-m', m, flg, ierr))
145:   n = 7
146:   PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', n, flg, ierr))
147:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))

149: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: !      Compute the matrix and right-hand-side vector that define
151: !      the linear system, Ax = b.
152: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

154: !  Create parallel matrix, specifying only its global dimensions.
155: !  When using MatCreate(), the matrix format can be specified at
156: !  runtime. Also, the parallel partitioning of the matrix is
157: !  determined by PETSc at runtime.

159:   PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
160:   PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m*n, m*n, ierr))
161:   PetscCallA(MatSetType(A, MATAIJ, ierr))
162:   PetscCallA(MatSetFromOptions(A, ierr))
163:   PetscCallA(MatMPIAIJSetPreallocation(A, 5_PETSC_INT_KIND, PETSC_NULL_INTEGER_ARRAY, 5_PETSC_INT_KIND, PETSC_NULL_INTEGER_ARRAY, ierr))
164:   PetscCallA(MatSeqAIJSetPreallocation(A, 5_PETSC_INT_KIND, PETSC_NULL_INTEGER_ARRAY, ierr))

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

170:   PetscCallA(MatGetOwnershipRange(A, Istart, Iend, ierr))

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

180:   do II = Istart, Iend - 1
181:     v = -1.0
182:     i = II/n
183:     j = II - i*n
184:     if (i > 0) then
185:       JJ = II - n
186:       PetscCallA(MatSetValues(A, 1_PETSC_INT_KIND, [II], 1_PETSC_INT_KIND, [JJ], [v], ADD_VALUES, ierr))
187:     end if
188:     if (i < m - 1) then
189:       JJ = II + n
190:       PetscCallA(MatSetValues(A, 1_PETSC_INT_KIND, [II], 1_PETSC_INT_KIND, [JJ], [v], ADD_VALUES, ierr))
191:     end if
192:     if (j > 0) then
193:       JJ = II - 1
194:       PetscCallA(MatSetValues(A, 1_PETSC_INT_KIND, [II], 1_PETSC_INT_KIND, [JJ], [v], ADD_VALUES, ierr))
195:     end if
196:     if (j < n - 1) then
197:       JJ = II + 1
198:       PetscCallA(MatSetValues(A, 1_PETSC_INT_KIND, [II], 1_PETSC_INT_KIND, [JJ], [v], ADD_VALUES, ierr))
199:     end if
200:     v = 4.0
201:     PetscCallA(MatSetValues(A, 1_PETSC_INT_KIND, [II], 1_PETSC_INT_KIND, [II], [v], ADD_VALUES, ierr))
202:   end do

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

209:   PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
210:   PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))

212: !  Create parallel vectors.
213: !   - Here, the parallel partitioning of the vector is determined by
214: !     PETSc at runtime.  We could also specify the local dimensions
215: !     if desired -- or use the more general routine VecCreate().
216: !   - When solving a linear system, the vectors and matrices MUST
217: !     be partitioned accordingly.  PETSc automatically generates
218: !     appropriately partitioned matrices and vectors when MatCreate()
219: !     and VecCreate() are used with the same communicator.
220: !   - Note: We form 1 vector from scratch and then duplicate as needed.

222:   PetscCallA(VecCreateFromOptions(PETSC_COMM_WORLD, PETSC_NULL_CHARACTER, 1_PETSC_INT_KIND, PETSC_DECIDE, m*n, u, ierr))
223:   PetscCallA(VecDuplicate(u, b, ierr))
224:   PetscCallA(VecDuplicate(b, x, ierr))

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

228:   PetscCallA(VecSet(u, one, ierr))
229:   PetscCallA(MatMult(A, u, b, ierr))

231: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232: !         Create the linear solver and set various options
233: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

235: !  Create linear solver context

237:   PetscCallA(KSPCreate(PETSC_COMM_WORLD, ksp, ierr))

239: !  Set operators. Here the matrix that defines the linear system
240: !  also serves as the matrix from which the preconditioner is constructed.

242:   PetscCallA(KSPSetOperators(ksp, A, A, ierr))

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

249:   PetscCallA(KSPGetPC(ksp, pc, ierr))
250:   PetscCallA(KSPSetTolerances(ksp, tol, PETSC_CURRENT_REAL, PETSC_CURRENT_REAL, PETSC_CURRENT_INTEGER, ierr))

252: !
253: !  Set a user-defined shell preconditioner if desired
254: !
255:   PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-user_defined_pc', user_defined_pc, ierr))

257:   if (user_defined_pc) then

259: !  (Required) Indicate to PETSc that we are using a shell preconditioner
260:     PetscCallA(PCSetType(pc, PCSHELL, ierr))

262: !  (Required) Set the user-defined routine for applying the preconditioner
263:     PetscCallA(PCShellSetApply(pc, SampleShellPCApply, ierr))

265: !  (Optional) Do any setup required for the preconditioner
266:     PetscCallA(PCShellSetSetUp(pc, SampleShellPCSetUp, ierr))

268: !  (Optional) Frees any objects we created for the preconditioner
269:     PetscCallA(PCShellSetDestroy(pc, SampleShellPCDestroy, ierr))

271:   else
272:     PetscCallA(PCSetType(pc, PCJACOBI, ierr))
273:   end if

275: !  Set runtime options, e.g.,
276: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
277: !  These options will override those specified above as long as
278: !  KSPSetFromOptions() is called _after_ any other customization
279: !  routines.

281:   PetscCallA(KSPSetFromOptions(ksp, ierr))

283: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
284: !                      Solve the linear system
285: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

287:   PetscCallA(KSPSolve(ksp, b, x, ierr))

289: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
290: !                     Check solution and clean up
291: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

293: !  Check the error

295:   PetscCallA(VecAXPY(x, neg_one, u, ierr))
296:   PetscCallA(VecNorm(x, NORM_2, norm, ierr))
297:   PetscCallA(KSPGetIterationNumber(ksp, its, ierr))

299:   if (rank == 0) then
300:     if (norm > 1.e-12) then
301:       write (6, 100) norm, its
302:     else
303:       write (6, 110) its
304:     end if
305:   end if
306: 100 format('Norm of error ', 1pe11.4, ' iterations ', i5)
307: 110 format('Norm of error < 1.e-12,iterations ', i5)

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

312:   PetscCallA(KSPDestroy(ksp, ierr))
313:   PetscCallA(VecDestroy(u, ierr))
314:   PetscCallA(VecDestroy(x, ierr))
315:   PetscCallA(VecDestroy(b, ierr))
316:   PetscCallA(MatDestroy(A, ierr))

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

320:   PetscCallA(PetscFinalize(ierr))
321: end program
322: !
323: !/*TEST
324: !
325: !   test:
326: !      nsize: 2
327: !      args: -ksp_view -user_defined_pc -ksp_gmres_cgs_refinement_type refine_always
328: !
329: !TEST*/