Actual source code: ex2f.F90

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: !
  2: !  Description: Solves a linear system in parallel with KSP (Fortran code).
  3: !               Also shows how to set a user-defined monitoring routine.
  4: !
  5: !
  6: !/*T
  7: !  Concepts: KSP^basic parallel example
  8: !  Concepts: KSP^setting a user-defined monitoring routine
  9: !  Processors: n
 10: !T*/
 11: !
 12: ! -----------------------------------------------------------------------

 14:       program main
 15:  #include <petsc/finclude/petscksp.h>
 16:       use petscksp
 17:       implicit none
 18: #if defined(PETSC_USING_F90) && !defined(PETSC_USE_FORTRANKIND)
 19:       external PETSC_NULL_FUNCTION
 20:       external KSPMONITORDEFAULT
 21:       external PETSCVIEWERANDFORMATDESTROY
 22: #endif
 23: !
 24: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 25: !                   Variable declarations
 26: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 27: !
 28: !  Variables:
 29: !     ksp     - linear solver context
 30: !     ksp      - Krylov subspace method context
 31: !     pc       - preconditioner context
 32: !     x, b, u  - approx solution, right-hand-side, exact solution vectors
 33: !     A        - matrix that defines linear system
 34: !     its      - iterations for convergence
 35: !     norm     - norm of error in solution
 36: !     rctx     - random number generator context
 37: !
 38: !  Note that vectors are declared as PETSc "Vec" objects.  These vectors
 39: !  are mathematical objects that contain more than just an array of
 40: !  double precision numbers. I.e., vectors in PETSc are not just
 41: !        double precision x(*).
 42: !  However, local vector data can be easily accessed via VecGetArray().
 43: !  See the Fortran section of the PETSc users manual for details.
 44: !
 45:       PetscReal  norm
 46:       PetscInt  i,j,II,JJ,m,n,its
 47:       PetscInt  Istart,Iend,ione
 48:       PetscErrorCode ierr
 49:       PetscMPIInt     rank,size
 50:       PetscBool   flg
 51:       PetscScalar v,one,neg_one
 52:       Vec         x,b,u
 53:       Mat         A
 54:       KSP         ksp
 55:       PetscRandom rctx
 56:       PetscViewerAndFormat vf;

 58: !  These variables are not currently used.
 59: !      PC          pc
 60: !      PCType      ptype
 61: !      PetscReal tol


 64: !  Note: Any user-defined Fortran routines (such as MyKSPMonitor)
 65: !  MUST be declared as external.

 67:       external MyKSPMonitor,MyKSPConverged

 69: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70: !                 Beginning of program
 71: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 73:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 74:       if (ierr .ne. 0) then
 75:         print*,'Unable to initialize PETSc'
 76:         stop
 77:       endif
 78:       m = 3
 79:       n = 3
 80:       one  = 1.0
 81:       neg_one = -1.0
 82:       ione    = 1
 83:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-m',m,flg,ierr)
 84:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
 85:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 86:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)

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

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

 98:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 99:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,ierr)
100:       call MatSetFromOptions(A,ierr)
101:       call MatSetUp(A,ierr)

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

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

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

117: !     Note: this uses the less common natural ordering that orders first
118: !     all the unknowns for x = h then for x = 2h etc; Hence you see JH = II +- n
119: !     instead of JJ = II +- m as you might expect. The more standard ordering
120: !     would first do all variables for y = h, then y = 2h etc.

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

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

151:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
152:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

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

164:       call VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,m*n,u,ierr)
165:       call VecSetFromOptions(u,ierr)
166:       call VecDuplicate(u,b,ierr)
167:       call VecDuplicate(b,x,ierr)

169: !  Set exact solution; then compute right-hand-side vector.
170: !  By default we use an exact solution of a vector with all
171: !  elements of 1.0;  Alternatively, using the runtime option
172: !  -random_sol forms a solution vector with random components.

174:       call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-random_exact_sol',flg,ierr)
175:       if (flg) then
176:          call PetscRandomCreate(PETSC_COMM_WORLD,rctx,ierr)
177:          call PetscRandomSetFromOptions(rctx,ierr)
178:          call VecSetRandom(u,rctx,ierr)
179:          call PetscRandomDestroy(rctx,ierr)
180:       else
181:          call VecSet(u,one,ierr)
182:       endif
183:       call MatMult(A,u,b,ierr)

185: !  View the exact solution vector if desired

187:       call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-view_exact_sol',flg,ierr)
188:       if (flg) then
189:          call VecView(u,PETSC_VIEWER_STDOUT_WORLD,ierr)
190:       endif

192: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193: !         Create the linear solver and set various options
194: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

196: !  Create linear solver context

198:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

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

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

205: !  Set linear solver defaults for this problem (optional).
206: !   - By extracting the KSP and PC contexts from the KSP context,
207: !     we can then directly directly call any KSP and PC routines
208: !     to set various options.
209: !   - The following four statements are optional; all of these
210: !     parameters could alternatively be specified at runtime via
211: !     KSPSetFromOptions(). All of these defaults can be
212: !     overridden at runtime, as indicated below.

214: !     We comment out this section of code since the Jacobi
215: !     preconditioner is not a good general default.

217: !      call KSPGetPC(ksp,pc,ierr)
218: !      ptype = PCJACOBI
219: !      call PCSetType(pc,ptype,ierr)
220: !      tol = 1.e-7
221: !      call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_REAL,PETSC_DEFAULT_REAL,PETSC_DEFAULT_INTEGER,ierr)

223: !  Set user-defined monitoring routine if desired

225:       call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_ksp_monitor',flg,ierr)
226:       if (flg) then
227:         call KSPMonitorSet(ksp,MyKSPMonitor,0,PETSC_NULL_FUNCTION,ierr)
228: !
229: !     Also use the default KSP monitor routine showing how it may be used from Fortran
230: !
231:         call PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,vf,ierr)
232:         call KSPMonitorSet(ksp,KSPMonitorDefault,vf,PetscViewerAndFormatDestroy,ierr)
233:       endif


236: !  Set runtime options, e.g.,
237: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
238: !  These options will override those specified above as long as
239: !  KSPSetFromOptions() is called _after_ any other customization
240: !  routines.

242:       call KSPSetFromOptions(ksp,ierr)

244: !  Set convergence test routine if desired

246:       call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_ksp_convergence',flg,ierr)
247:       if (flg) then
248:         call KSPSetConvergenceTest(ksp,MyKSPConverged,0,PETSC_NULL_FUNCTION,ierr)
249:       endif
250: !
251: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252: !                      Solve the linear system
253: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

257: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258: !                     Check solution and clean up
259: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

261: !  Check the error
262:       call VecAXPY(x,neg_one,u,ierr)
263:       call VecNorm(x,NORM_2,norm,ierr)
264:       call KSPGetIterationNumber(ksp,its,ierr)
265:       if (rank .eq. 0) then
266:         if (norm .gt. 1.e-12) then
267:            write(6,100) norm,its
268:         else
269:            write(6,110) its
270:         endif
271:       endif
272:   100 format('Norm of error ',e11.4,' iterations ',i5)
273:   110 format('Norm of error < 1.e-12 iterations ',i5)

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

278:       call KSPDestroy(ksp,ierr)
279:       call VecDestroy(u,ierr)
280:       call VecDestroy(x,ierr)
281:       call VecDestroy(b,ierr)
282:       call MatDestroy(A,ierr)

284: !  Always call PetscFinalize() before exiting a program.  This routine
285: !    - finalizes the PETSc libraries as well as MPI
286: !    - provides summary and diagnostic information if certain runtime
287: !      options are chosen (e.g., -log_view).  See PetscFinalize()
288: !      manpage for more information.

290:       call PetscFinalize(ierr)
291:       end

293: ! --------------------------------------------------------------
294: !
295: !  MyKSPMonitor - This is a user-defined routine for monitoring
296: !  the KSP iterative solvers.
297: !
298: !  Input Parameters:
299: !    ksp   - iterative context
300: !    n     - iteration number
301: !    rnorm - 2-norm (preconditioned) residual value (may be estimated)
302: !    dummy - optional user-defined monitor context (unused here)
303: !
304:       subroutine MyKSPMonitor(ksp,n,rnorm,dummy,ierr)
305:       use petscksp
306:       implicit none

308:       KSP              ksp
309:       Vec              x
310:       PetscErrorCode ierr
311:       PetscInt n,dummy
312:       PetscMPIInt rank
313:       PetscReal rnorm

315: !  Build the solution vector
316:       call KSPBuildSolution(ksp,PETSC_NULL_VEC,x,ierr)

318: !  Write the solution vector and residual norm to stdout
319: !   - Note that the parallel viewer PETSC_VIEWER_STDOUT_WORLD
320: !     handles data from multiple processors so that the
321: !     output is not jumbled.

323:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
324:       if (rank .eq. 0) write(6,100) n
325:       call VecView(x,PETSC_VIEWER_STDOUT_WORLD,ierr)
326:       if (rank .eq. 0) write(6,200) n,rnorm

328:  100  format('iteration ',i5,' solution vector:')
329:  200  format('iteration ',i5,' residual norm ',e11.4)
330:       0
331:       end

333: ! --------------------------------------------------------------
334: !
335: !  MyKSPConverged - This is a user-defined routine for testing
336: !  convergence of the KSP iterative solvers.
337: !
338: !  Input Parameters:
339: !    ksp   - iterative context
340: !    n     - iteration number
341: !    rnorm - 2-norm (preconditioned) residual value (may be estimated)
342: !    dummy - optional user-defined monitor context (unused here)
343: !
344:       subroutine MyKSPConverged(ksp,n,rnorm,flag,dummy,ierr)
345:       use petscksp
346:       implicit none

348:       KSP              ksp
349:       PetscErrorCode ierr
350:       PetscInt n,dummy
351:       KSPConvergedReason flag
352:       PetscReal rnorm

354:       if (rnorm .le. .05) then
355:         flag = 1
356:       else
357:         flag = 0
358:       endif
359:       0

361:       end

363: !/*TEST
364: !
365: !   test:
366: !      nsize: 2
367: !      args: -pc_type jacobi -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
368: !
369: !   test:
370: !      suffix: 2
371: !      nsize: 2
372: !      args: -pc_type jacobi -my_ksp_monitor -ksp_gmres_cgs_refinement_type refine_always
373: !
374: !TEST*/