Actual source code: ex1f.F

slepc-3.20.1 2023-11-27
Report Typos and Errors
  1: !
  2: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3: !  SLEPc - Scalable Library for Eigenvalue Problem Computations
  4: !  Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
  5: !
  6: !  This file is part of SLEPc.
  7: !  SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: !
 10: !  Program usage: mpiexec -n <np> ./ex1f [-help] [-n <n>] [all SLEPc options]
 11: !
 12: !  Description: Simple example that solves an eigensystem with the EPS object.
 13: !  The standard symmetric eigenvalue problem to be solved corresponds to the
 14: !  Laplacian operator in 1 dimension.
 15: !
 16: !  The command line options are:
 17: !    -n <n>, where <n> = number of grid points = matrix size
 18: !
 19: !  Notes:
 20: !  - The free-form equivalent of this example is ex1f90.F90.
 21: !  - Does not check errors, see ex1f90.F90 for error checking (recommended).
 22: ! ----------------------------------------------------------------------
 23: !
 24:       program main
 25: #include <slepc/finclude/slepceps.h>
 26:       use slepceps
 27:       implicit none

 29: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 30: !     Declarations
 31: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32: !
 33: !  Variables:
 34: !     A     operator matrix
 35: !     eps   eigenproblem solver context

 37:       Mat            A
 38:       EPS            eps
 39:       EPSType        tname
 40:       PetscReal      tol, error
 41:       PetscScalar    kr, ki
 42:       Vec            xr, xi
 43:       PetscInt       n, i, Istart, Iend
 44:       PetscInt       nev, maxit, its, nconv
 45:       PetscInt       col(3)
 46:       PetscInt       i1,i2,i3
 47:       PetscMPIInt    rank
 48:       PetscErrorCode ierr
 49:       PetscBool      flg
 50:       PetscScalar    value(3)

 52: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 53: !     Beginning of program
 54: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 56:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 57:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 58:       n = 30
 59:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,  &
 60:      &                        '-n',n,flg,ierr)

 62:       if (rank .eq. 0) then
 63:         write(*,100) n
 64:       endif
 65:  100  format (/'1-D Laplacian Eigenproblem, n =',I3,' (Fortran)')

 67: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68: !     Compute the operator matrix that defines the eigensystem, Ax=kx
 69: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 71:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 72:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 73:       call MatSetFromOptions(A,ierr)
 74:       call MatSetUp(A,ierr)

 76:       i1 = 1
 77:       i2 = 2
 78:       i3 = 3
 79:       call MatGetOwnershipRange(A,Istart,Iend,ierr)
 80:       if (Istart .eq. 0) then
 81:         i = 0
 82:         col(1) = 0
 83:         col(2) = 1
 84:         value(1) =  2.0
 85:         value(2) = -1.0
 86:         call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
 87:         Istart = Istart+1
 88:       endif
 89:       if (Iend .eq. n) then
 90:         i = n-1
 91:         col(1) = n-2
 92:         col(2) = n-1
 93:         value(1) = -1.0
 94:         value(2) =  2.0
 95:         call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
 96:         Iend = Iend-1
 97:       endif
 98:       value(1) = -1.0
 99:       value(2) =  2.0
100:       value(3) = -1.0
101:       do i=Istart,Iend-1
102:         col(1) = i-1
103:         col(2) = i
104:         col(3) = i+1
105:         call MatSetValues(A,i1,i,i3,col,value,INSERT_VALUES,ierr)
106:       enddo

108:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
109:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

111:       call MatCreateVecs(A,xr,xi,ierr)

113: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: !     Create the eigensolver and display info
115: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

117: !     ** Create eigensolver context
118:       call EPSCreate(PETSC_COMM_WORLD,eps,ierr)

120: !     ** Set operators. In this case, it is a standard eigenvalue problem
121:       call EPSSetOperators(eps,A,PETSC_NULL_MAT,ierr)
122:       call EPSSetProblemType(eps,EPS_HEP,ierr)

124: !     ** Set solver parameters at runtime
125:       call EPSSetFromOptions(eps,ierr)

127: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: !     Solve the eigensystem
129: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

131:       call EPSSolve(eps,ierr)
132:       call EPSGetIterationNumber(eps,its,ierr)
133:       if (rank .eq. 0) then
134:         write(*,110) its
135:       endif
136:  110  format (/' Number of iterations of the method:',I4)

138: !     ** Optional: Get some information from the solver and display it
139:       call EPSGetType(eps,tname,ierr)
140:       if (rank .eq. 0) then
141:         write(*,120) tname
142:       endif
143:  120  format (' Solution method: ',A)
144:       call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,                 &
145:      &                      PETSC_NULL_INTEGER,ierr)
146:       if (rank .eq. 0) then
147:         write(*,130) nev
148:       endif
149:  130  format (' Number of requested eigenvalues:',I2)
150:       call EPSGetTolerances(eps,tol,maxit,ierr)
151:       if (rank .eq. 0) then
152:         write(*,140) tol, maxit
153:       endif
154:  140  format (' Stopping condition: tol=',1P,E11.4,', maxit=',I4)

156: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: !     Display solution and clean up
158: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

160: !     ** Get number of converged eigenpairs
161:       call EPSGetConverged(eps,nconv,ierr)
162:       if (rank .eq. 0) then
163:         write(*,150) nconv
164:       endif
165:  150  format (' Number of converged eigenpairs:',I2/)

167: !     ** Display eigenvalues and relative errors
168:       if (nconv.gt.0) then
169:         if (rank .eq. 0) then
170:           write(*,*) '         k          ||Ax-kx||/||kx||'
171:           write(*,*) ' ----------------- ------------------'
172:         endif
173:         do i=0,nconv-1
174: !         ** Get converged eigenpairs: i-th eigenvalue is stored in kr
175: !         ** (real part) and ki (imaginary part)
176:           call EPSGetEigenpair(eps,i,kr,ki,xr,xi,ierr)

178: !         ** Compute the relative error associated to each eigenpair
179:           call EPSComputeError(eps,i,EPS_ERROR_RELATIVE,error,ierr)
180:           if (rank .eq. 0) then
181:             write(*,160) PetscRealPart(kr), error
182:           endif
183:  160      format (1P,'   ',E12.4,'       ',E12.4)

185:         enddo
186:         if (rank .eq. 0) then
187:           write(*,*)
188:         endif
189:       endif

191: !     ** Free work space
192:       call EPSDestroy(eps,ierr)
193:       call MatDestroy(A,ierr)
194:       call VecDestroy(xr,ierr)
195:       call VecDestroy(xi,ierr)

197:       call SlepcFinalize(ierr)
198:       end