Actual source code: ex1f.F

  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:       implicit none

 13: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 14: !                    Include files
 15: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 16: !
 17: !  This program uses CPP for preprocessing, as indicated by the use of
 18: !  PETSc include files in the directory petsc/include/finclude.  This
 19: !  convention enables use of the CPP preprocessor, which allows the use
 20: !  of the #include statements that define PETSc objects and variables.
 21: !
 22: !  Use of the conventional Fortran include statements is also supported
 23: !  In this case, the PETsc include files are located in the directory
 24: !  petsc/include/foldinclude.
 25: !
 26: !  Since one must be very careful to include each file no more than once
 27: !  in a Fortran routine, application programmers must exlicitly list
 28: !  each file needed for the various PETSc components within their
 29: !  program (unlike the C/C++ interface).
 30: !
 31: !  See the Fortran section of the PETSc users manual for details.
 32: !
 33: !  The following include statements are required for KSP Fortran programs:
 34: !     petsc.h       - base PETSc routines
 35: !     petscvec.h    - vectors
 36: !     petscmat.h    - matrices
 37: !     petscksp.h    - Krylov subspace methods
 38: !     petscpc.h     - preconditioners
 39: !  Other include statements may be needed if using additional PETSc
 40: !  routines in a Fortran program, e.g.,
 41: !     petscviewer.h - viewers
 42: !     petscis.h     - index sets
 43: !
 44:  #include include/finclude/petsc.h
 45:  #include include/finclude/petscvec.h
 46:  #include include/finclude/petscmat.h
 47:  #include include/finclude/petscksp.h
 48:  #include include/finclude/petscpc.h
 49: !
 50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51: !                   Variable declarations
 52: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 53: !
 54: !  Variables:
 55: !     ksp     - linear solver context
 56: !     ksp      - Krylov subspace method context
 57: !     pc       - preconditioner context
 58: !     x, b, u  - approx solution, right-hand-side, exact solution vectors
 59: !     A        - matrix that defines linear system
 60: !     its      - iterations for convergence
 61: !     norm     - norm of error in solution
 62: !
 63:       Vec              x,b,u
 64:       Mat              A
 65:       KSP              ksp
 66:       PC               pc
 67:       double precision norm,tol
 68:       PetscErrorCode   ierr
 69:       PetscInt i,n,col(3),its,i1,i2,i3
 70:       PetscTruth flg
 71:       PetscMPIInt size,rank
 72:       PetscScalar      none,one,value(3)

 74: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75: !                 Beginning of program
 76: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 78:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 79:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
 80:       if (size .ne. 1) then
 81:          call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 82:          if (rank .eq. 0) then
 83:             write(6,*) 'This is a uniprocessor example only!'
 84:          endif
 85:          SETERRQ(1,' ',ierr)
 86:       endif
 87:       none = -1.0
 88:       one  = 1.0
 89:       n    = 10
 90:       i1 = 1
 91:       i2 = 2
 92:       i3 = 3
 93:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)

 95: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96: !         Compute the matrix and right-hand-side vector that define
 97: !         the linear system, Ax = b.
 98: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

103:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
104:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
105:       call MatSetFromOptions(A,ierr)

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

111:       value(1) = -1.0
112:       value(2) = 2.0
113:       value(3) = -1.0
114:       do 50 i=1,n-2
115:          col(1) = i-1
116:          col(2) = i
117:          col(3) = i+1
118:          call MatSetValues(A,i1,i,i3,col,value,INSERT_VALUES,ierr)
119:   50  continue
120:       i = n - 1
121:       col(1) = n - 2
122:       col(2) = n - 1
123:       call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
124:       i = 0
125:       col(1) = 0
126:       col(2) = 1
127:       value(1) = 2.0
128:       value(2) = -1.0
129:       call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
130:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
131:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

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

136:       call VecCreate(PETSC_COMM_WORLD,x,ierr)
137:       call VecSetSizes(x,PETSC_DECIDE,n,ierr)
138:       call VecSetFromOptions(x,ierr)
139:       call VecDuplicate(x,b,ierr)
140:       call VecDuplicate(x,u,ierr)

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

144:       call VecSet(u,one,ierr)
145:       call MatMult(A,u,b,ierr)

147: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: !          Create the linear solver and set various options
149: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

151: !  Create linear solver context

153:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

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

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

160: !  Set linear solver defaults for this problem (optional).
161: !   - By extracting the KSP and PC contexts from the KSP context,
162: !     we can then directly directly call any KSP and PC routines
163: !     to set various options.
164: !   - The following four statements are optional; all of these
165: !     parameters could alternatively be specified at runtime via
166: !     KSPSetFromOptions();

168:       call KSPGetPC(ksp,pc,ierr)
169:       call PCSetType(pc,PCJACOBI,ierr)
170:       tol = 1.d-7
171:       call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_DOUBLE_PRECISION,     &
172:      &     PETSC_DEFAULT_DOUBLE_PRECISION,PETSC_DEFAULT_INTEGER,ierr)

174: !  Set runtime options, e.g.,
175: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
176: !  These options will override those specified above as long as
177: !  KSPSetFromOptions() is called _after_ any other customization
178: !  routines.

180:       call KSPSetFromOptions(ksp,ierr)

182: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183: !                      Solve the linear system
184: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

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

190:       call KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD,ierr)

192: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193: !                      Check solution and clean up
194: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

196: !  Check the error

198:       call VecAXPY(x,none,u,ierr)
199:       call VecNorm(x,NORM_2,norm,ierr)
200:       call KSPGetIterationNumber(ksp,its,ierr)
201:       if (norm .gt. 1.e-12) then
202:         write(6,100) norm,its
203:       else
204:         write(6,200) its
205:       endif
206:  100  format('Norm of error = ',e10.4,',  Iterations = ',i5)
207:  200  format('Norm of error < 1.e-12,Iterations = ',i5)

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

212:       call VecDestroy(x,ierr)
213:       call VecDestroy(u,ierr)
214:       call VecDestroy(b,ierr)
215:       call MatDestroy(A,ierr)
216:       call KSPDestroy(ksp,ierr)
217:       call PetscFinalize(ierr)

219:       end