Actual source code: ex1f.F

  1: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  2: !  SLEPc - Scalable Library for Eigenvalue Problem Computations
  3: !  Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
  4: !
  5: !  This file is part of SLEPc.
  6: !
  7: !  SLEPc is free software: you can redistribute it and/or modify it under  the
  8: !  terms of version 3 of the GNU Lesser General Public License as published by
  9: !  the Free Software Foundation.
 10: !
 11: !  SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 12: !  WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 13: !  FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 14: !  more details.
 15: !
 16: !  You  should have received a copy of the GNU Lesser General  Public  License
 17: !  along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 18: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 19: !
 20: !  Program usage: mpirun -np n ex1f [-help] [-n <n>] [all SLEPc options]
 21: !
 22: !  Description: Simple example that solves an eigensystem with the EPS object.
 23: !  The standard symmetric eigenvalue problem to be solved corresponds to the
 24: !  Laplacian operator in 1 dimension.
 25: !
 26: !  The command line options are:
 27: !    -n <n>, where <n> = number of grid points = matrix size
 28: !
 29: ! ----------------------------------------------------------------------
 30: !
 31:       program main
 32:       implicit none

 34: #include "finclude/petscsys.h"
 35: #include "finclude/petscvec.h"
 36: #include "finclude/petscmat.h"
 37:  #include finclude/slepcsys.h
 38:  #include finclude/slepceps.h

 40: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41: !     Declarations
 42: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 43: !
 44: !  Variables:
 45: !     A     operator matrix
 46: !     eps   eigenproblem solver context

 48:       Mat            A
 49:       EPS            eps
 50:       EPSType        tname
 51:       PetscReal      tol, error
 52:       PetscScalar    kr, ki
 53:       PetscInt       n, i, Istart, Iend
 54:       PetscInt       nev, maxit, its, nconv
 55:       PetscInt       col(3)
 56:       PetscInt       i1,i2,i3
 57:       PetscMPIInt    rank
 58:       PetscErrorCode ierr
 59:       PetscTruth     flg
 60:       PetscScalar    value(3)

 62: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63: !     Beginning of program
 64: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 66:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 67:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 68:       n = 30
 69:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)

 71:       if (rank .eq. 0) then
 72:         write(*,100) n
 73:       endif
 74:  100  format (/'1-D Laplacian Eigenproblem, n =',I3,' (Fortran)')

 76: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 77: !     Compute the operator matrix that defines the eigensystem, Ax=kx
 78: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 80:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 81:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 82:       call MatSetFromOptions(A,ierr)

 84:       i1 = 1
 85:       i2 = 2
 86:       i3 = 3
 87:       call MatGetOwnershipRange(A,Istart,Iend,ierr)
 88:       if (Istart .eq. 0) then
 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:         Istart = Istart+1
 96:       endif
 97:       if (Iend .eq. n) then
 98:         i = n-1
 99:         col(1) = n-2
100:         col(2) = n-1
101:         value(1) = -1.0
102:         value(2) =  2.0
103:         call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
104:         Iend = Iend-1
105:       endif
106:       value(1) = -1.0
107:       value(2) =  2.0
108:       value(3) = -1.0
109:       do i=Istart,Iend-1
110:         col(1) = i-1
111:         col(2) = i
112:         col(3) = i+1
113:         call MatSetValues(A,i1,i,i3,col,value,INSERT_VALUES,ierr)
114:       enddo

116:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
117:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

119: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: !     Create the eigensolver and display info
121: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

123: !     ** Create eigensolver context
124:       call EPSCreate(PETSC_COMM_WORLD,eps,ierr)

126: !     ** Set operators. In this case, it is a standard eigenvalue problem
127:       call EPSSetOperators(eps,A,PETSC_NULL_OBJECT,ierr)
128:       call EPSSetProblemType(eps,EPS_HEP,ierr)

130: !     ** Set solver parameters at runtime
131:       call EPSSetFromOptions(eps,ierr)

133: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: !     Solve the eigensystem
135: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

162: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163: !     Display solution and clean up
164: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

166: !     ** Get number of converged eigenpairs
167:       call EPSGetConverged(eps,nconv,ierr)
168:       if (rank .eq. 0) then
169:         write(*,150) nconv
170:       endif
171:  150  format (' Number of converged eigenpairs:',I2/)

173: !     ** Display eigenvalues and relative errors
174:       if (nconv.gt.0) then
175:         if (rank .eq. 0) then
176:           write(*,*) '         k          ||Ax-kx||/||kx||'
177:           write(*,*) ' ----------------- ------------------'
178:         endif
179:         do i=0,nconv-1
180: !         ** Get converged eigenpairs: i-th eigenvalue is stored in kr
181: !         ** (real part) and ki (imaginary part)
182:           call EPSGetEigenpair(eps,i,kr,ki,PETSC_NULL_OBJECT,
183:      +                         PETSC_NULL_OBJECT,ierr)

185: !         ** Compute the relative error associated to each eigenpair
186:           call EPSComputeRelativeError(eps,i,error,ierr)
187:           if (rank .eq. 0) then
188:             write(*,160) PetscRealPart(kr), error
189:           endif
190:  160      format (1P,'   ',E12.4,'       ',E12.4)

192:         enddo
193:         if (rank .eq. 0) then
194:           write(*,*)
195:         endif
196:       endif

198: !     ** Free work space
199:       call EPSDestroy(eps,ierr)
200:       call MatDestroy(A,ierr)

202:       call SlepcFinalize(ierr)
203:       end