Actual source code: ex1f.F90

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: !
  2: !   Description: Solves a tridiagonal linear system with KSP.
  3: !
  4: !/*T
  5: !   Concepts: KSP^solving a system of linear equations
  6: !   Processors: 1
  7: !T*/
  8: ! -----------------------------------------------------------------------

 10:       program main
 11:  #include <petsc/finclude/petscksp.h>
 12:       use petscksp
 13:       implicit none

 15: !
 16: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 17: !                   Variable declarations
 18: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 19: !
 20: !  Variables:
 21: !     ksp     - linear solver context
 22: !     ksp      - Krylov subspace method context
 23: !     pc       - preconditioner context
 24: !     x, b, u  - approx solution, right-hand-side, exact solution vectors
 25: !     A        - matrix that defines linear system
 26: !     its      - iterations for convergence
 27: !     norm     - norm of error in solution
 28: !
 29:       Vec              x,b,u
 30:       Mat              A
 31:       KSP              ksp
 32:       PC               pc
 33:       PetscReal        norm,tol
 34:       PetscErrorCode   ierr
 35:       PetscInt i,n,col(3),its,i1,i2,i3
 36:       PetscBool  flg
 37:       PetscMPIInt size
 38:       PetscScalar      none,one,value(3)

 40: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41: !                 Beginning of program
 42: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 44:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 45:       if (ierr .ne. 0) then
 46:         print*,'Unable to initialize PETSc'
 47:         stop
 48:       endif
 49:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
 50:       if (size .ne. 1) then; SETERRA(PETSC_COMM_WORLD,1,'This is a uniprocessor example only'); endif
 51:       none = -1.0
 52:       one  = 1.0
 53:       n    = 10
 54:       i1 = 1
 55:       i2 = 2
 56:       i3 = 3
 57:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr)

 59: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60: !         Compute the matrix and right-hand-side vector that define
 61: !         the linear system, Ax = b.
 62: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 64: !  Create matrix.  When using MatCreate(), the matrix format can
 65: !  be specified at runtime.

 67:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 68:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 69:       call MatSetFromOptions(A,ierr)
 70:       call MatSetUp(A,ierr)

 72: !  Assemble matrix.
 73: !   - Note that MatSetValues() uses 0-based row and column numbers
 74: !     in Fortran as well as in C (as set here in the array "col").

 76:       value(1) = -1.0
 77:       value(2) = 2.0
 78:       value(3) = -1.0
 79:       do 50 i=1,n-2
 80:          col(1) = i-1
 81:          col(2) = i
 82:          col(3) = i+1
 83:          call MatSetValues(A,i1,i,i3,col,value,INSERT_VALUES,ierr)
 84:   50  continue
 85:       i = n - 1
 86:       col(1) = n - 2
 87:       col(2) = n - 1
 88:       call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
 89:       i = 0
 90:       col(1) = 0
 91:       col(2) = 1
 92:       value(1) = 2.0
 93:       value(2) = -1.0
 94:       call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
 95:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
 96:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

 98: !  Create vectors.  Note that we form 1 vector from scratch and
 99: !  then duplicate as needed.

101:       call VecCreate(PETSC_COMM_WORLD,x,ierr)
102:       call VecSetSizes(x,PETSC_DECIDE,n,ierr)
103:       call VecSetFromOptions(x,ierr)
104:       call VecDuplicate(x,b,ierr)
105:       call VecDuplicate(x,u,ierr)

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

109:       call VecSet(u,one,ierr)
110:       call MatMult(A,u,b,ierr)

112: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: !          Create the linear solver and set various options
114: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

116: !  Create linear solver context

118:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

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

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

125: !  Set linear solver defaults for this problem (optional).
126: !   - By extracting the KSP and PC contexts from the KSP context,
127: !     we can then directly directly call any KSP and PC routines
128: !     to set various options.
129: !   - The following four statements are optional; all of these
130: !     parameters could alternatively be specified at runtime via
131: !     KSPSetFromOptions();

133:       call KSPGetPC(ksp,pc,ierr)
134:       call PCSetType(pc,PCJACOBI,ierr)
135:       tol = .0000001
136:       call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_REAL,                         &
137:      &     PETSC_DEFAULT_REAL,PETSC_DEFAULT_INTEGER,ierr)

139: !  Set runtime options, e.g.,
140: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
141: !  These options will override those specified above as long as
142: !  KSPSetFromOptions() is called _after_ any other customization
143: !  routines.

145:       call KSPSetFromOptions(ksp,ierr)

147: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: !                      Solve the linear system
149: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

153: !  View solver info; we could instead use the option -ksp_view

155:       call KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD,ierr)

157: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: !                      Check solution and clean up
159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

161: !  Check the error

163:       call VecAXPY(x,none,u,ierr)
164:       call VecNorm(x,NORM_2,norm,ierr)
165:       call KSPGetIterationNumber(ksp,its,ierr)
166:       if (norm .gt. 1.e-12) then
167:         write(6,100) norm,its
168:       else
169:         write(6,200) its
170:       endif
171:  100  format('Norm of error ',e11.4,',  Iterations = ',i5)
172:  200  format('Norm of error < 1.e-12, Iterations = ',i5)

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

177:       call VecDestroy(x,ierr)
178:       call VecDestroy(u,ierr)
179:       call VecDestroy(b,ierr)
180:       call MatDestroy(A,ierr)
181:       call KSPDestroy(ksp,ierr)
182:       call PetscFinalize(ierr)

184:       end

186: !/*TEST
187: !
188: !     test:
189: !       args: -ksp_monitor_short
190: !
191: !TEST*/