Actual source code: ex1f.F
slepc-3.7.0 2016-05-16
1: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2: ! SLEPc - Scalable Library for Eigenvalue Problem Computations
3: ! Copyright (c) 2002-2016, Universitat 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: mpiexec -n <np> ./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 <petsc/finclude/petscsys.h>
35: #include <petsc/finclude/petscvec.h>
36: #include <petsc/finclude/petscmat.h>
37: #include <slepc/finclude/slepcsys.h>
38: #include <slepc/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: Vec xr, xi
54: PetscInt n, i, Istart, Iend
55: PetscInt nev, maxit, its, nconv
56: PetscInt col(3)
57: PetscInt i1,i2,i3
58: PetscMPIInt rank
59: PetscErrorCode ierr
60: PetscBool flg
61: PetscScalar value(3)
63: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: ! Beginning of program
65: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
68: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
69: n = 30
70: call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER, &
71: & '-n',n,flg,ierr)
73: if (rank .eq. 0) then
74: write(*,100) n
75: endif
76: 100 format (/'1-D Laplacian Eigenproblem, n =',I3,' (Fortran)')
78: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: ! Compute the operator matrix that defines the eigensystem, Ax=kx
80: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: call MatCreate(PETSC_COMM_WORLD,A,ierr)
83: call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
84: call MatSetFromOptions(A,ierr)
85: call MatSetUp(A,ierr)
87: i1 = 1
88: i2 = 2
89: i3 = 3
90: call MatGetOwnershipRange(A,Istart,Iend,ierr)
91: if (Istart .eq. 0) then
92: i = 0
93: col(1) = 0
94: col(2) = 1
95: value(1) = 2.0
96: value(2) = -1.0
97: call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
98: Istart = Istart+1
99: endif
100: if (Iend .eq. n) then
101: i = n-1
102: col(1) = n-2
103: col(2) = n-1
104: value(1) = -1.0
105: value(2) = 2.0
106: call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
107: Iend = Iend-1
108: endif
109: value(1) = -1.0
110: value(2) = 2.0
111: value(3) = -1.0
112: do i=Istart,Iend-1
113: col(1) = i-1
114: col(2) = i
115: col(3) = i+1
116: call MatSetValues(A,i1,i,i3,col,value,INSERT_VALUES,ierr)
117: enddo
119: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
120: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
122: call MatCreateVecs(A,xr,xi,ierr)
124: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: ! Create the eigensolver and display info
126: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: ! ** Create eigensolver context
129: call EPSCreate(PETSC_COMM_WORLD,eps,ierr)
131: ! ** Set operators. In this case, it is a standard eigenvalue problem
132: call EPSSetOperators(eps,A,PETSC_NULL_OBJECT,ierr)
133: call EPSSetProblemType(eps,EPS_HEP,ierr)
135: ! ** Set solver parameters at runtime
136: call EPSSetFromOptions(eps,ierr)
138: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: ! Solve the eigensystem
140: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: call EPSSolve(eps,ierr)
143: call EPSGetIterationNumber(eps,its,ierr)
144: if (rank .eq. 0) then
145: write(*,110) its
146: endif
147: 110 format (/' Number of iterations of the method:',I4)
149: ! ** Optional: Get some information from the solver and display it
150: call EPSGetType(eps,tname,ierr)
151: if (rank .eq. 0) then
152: write(*,120) tname
153: endif
154: 120 format (' Solution method: ',A)
155: call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER, &
156: & PETSC_NULL_INTEGER,ierr)
157: if (rank .eq. 0) then
158: write(*,130) nev
159: endif
160: 130 format (' Number of requested eigenvalues:',I2)
161: call EPSGetTolerances(eps,tol,maxit,ierr)
162: if (rank .eq. 0) then
163: write(*,140) tol, maxit
164: endif
165: 140 format (' Stopping condition: tol=',1P,E10.4,', maxit=',I4)
167: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168: ! Display solution and clean up
169: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171: ! ** Get number of converged eigenpairs
172: call EPSGetConverged(eps,nconv,ierr)
173: if (rank .eq. 0) then
174: write(*,150) nconv
175: endif
176: 150 format (' Number of converged eigenpairs:',I2/)
178: ! ** Display eigenvalues and relative errors
179: if (nconv.gt.0) then
180: if (rank .eq. 0) then
181: write(*,*) ' k ||Ax-kx||/||kx||'
182: write(*,*) ' ----------------- ------------------'
183: endif
184: do i=0,nconv-1
185: ! ** Get converged eigenpairs: i-th eigenvalue is stored in kr
186: ! ** (real part) and ki (imaginary part)
187: call EPSGetEigenpair(eps,i,kr,ki,xr,xi,ierr)
189: ! ** Compute the relative error associated to each eigenpair
190: call EPSComputeError(eps,i,EPS_ERROR_RELATIVE,error,ierr)
191: if (rank .eq. 0) then
192: write(*,160) PetscRealPart(kr), error
193: endif
194: 160 format (1P,' ',E12.4,' ',E12.4)
196: enddo
197: if (rank .eq. 0) then
198: write(*,*)
199: endif
200: endif
202: ! ** Free work space
203: call EPSDestroy(eps,ierr)
204: call MatDestroy(A,ierr)
205: call VecDestroy(xr,ierr)
206: call VecDestroy(xi,ierr)
208: call SlepcFinalize(ierr)
209: end