Actual source code: ex1f.F
slepc-3.6.1 2015-09-03
1: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2: ! SLEPc - Scalable Library for Eigenvalue Problem Computations
3: ! Copyright (c) 2002-2015, 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: 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 <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_CHARACTER,'-n',n,flg,ierr)
72: if (rank .eq. 0) then
73: write(*,100) n
74: endif
75: 100 format (/'1-D Laplacian Eigenproblem, n =',I3,' (Fortran)')
77: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: ! Compute the operator matrix that defines the eigensystem, Ax=kx
79: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: call MatCreate(PETSC_COMM_WORLD,A,ierr)
82: call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
83: call MatSetFromOptions(A,ierr)
84: call MatSetUp(A,ierr)
86: i1 = 1
87: i2 = 2
88: i3 = 3
89: call MatGetOwnershipRange(A,Istart,Iend,ierr)
90: if (Istart .eq. 0) then
91: i = 0
92: col(1) = 0
93: col(2) = 1
94: value(1) = 2.0
95: value(2) = -1.0
96: call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
97: Istart = Istart+1
98: endif
99: if (Iend .eq. n) then
100: i = n-1
101: col(1) = n-2
102: col(2) = n-1
103: value(1) = -1.0
104: value(2) = 2.0
105: call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
106: Iend = Iend-1
107: endif
108: value(1) = -1.0
109: value(2) = 2.0
110: value(3) = -1.0
111: do i=Istart,Iend-1
112: col(1) = i-1
113: col(2) = i
114: col(3) = i+1
115: call MatSetValues(A,i1,i,i3,col,value,INSERT_VALUES,ierr)
116: enddo
118: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
119: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
121: call MatCreateVecs(A,xr,xi,ierr)
123: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: ! Create the eigensolver and display info
125: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: ! ** Create eigensolver context
128: call EPSCreate(PETSC_COMM_WORLD,eps,ierr)
130: ! ** Set operators. In this case, it is a standard eigenvalue problem
131: call EPSSetOperators(eps,A,PETSC_NULL_OBJECT,ierr)
132: call EPSSetProblemType(eps,EPS_HEP,ierr)
134: ! ** Set solver parameters at runtime
135: call EPSSetFromOptions(eps,ierr)
137: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: ! Solve the eigensystem
139: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: call EPSSolve(eps,ierr)
142: call EPSGetIterationNumber(eps,its,ierr)
143: if (rank .eq. 0) then
144: write(*,110) its
145: endif
146: 110 format (/' Number of iterations of the method:',I4)
148: ! ** Optional: Get some information from the solver and display it
149: call EPSGetType(eps,tname,ierr)
150: if (rank .eq. 0) then
151: write(*,120) tname
152: endif
153: 120 format (' Solution method: ',A)
154: call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER, &
155: & PETSC_NULL_INTEGER,ierr)
156: if (rank .eq. 0) then
157: write(*,130) nev
158: endif
159: 130 format (' Number of requested eigenvalues:',I2)
160: call EPSGetTolerances(eps,tol,maxit,ierr)
161: if (rank .eq. 0) then
162: write(*,140) tol, maxit
163: endif
164: 140 format (' Stopping condition: tol=',1P,E10.4,', maxit=',I4)
166: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: ! Display solution and clean up
168: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170: ! ** Get number of converged eigenpairs
171: call EPSGetConverged(eps,nconv,ierr)
172: if (rank .eq. 0) then
173: write(*,150) nconv
174: endif
175: 150 format (' Number of converged eigenpairs:',I2/)
177: ! ** Display eigenvalues and relative errors
178: if (nconv.gt.0) then
179: if (rank .eq. 0) then
180: write(*,*) ' k ||Ax-kx||/||kx||'
181: write(*,*) ' ----------------- ------------------'
182: endif
183: do i=0,nconv-1
184: ! ** Get converged eigenpairs: i-th eigenvalue is stored in kr
185: ! ** (real part) and ki (imaginary part)
186: call EPSGetEigenpair(eps,i,kr,ki,xr,xi,ierr)
188: ! ** Compute the relative error associated to each eigenpair
189: call EPSComputeError(eps,i,EPS_ERROR_RELATIVE,error,ierr)
190: if (rank .eq. 0) then
191: write(*,160) PetscRealPart(kr), error
192: endif
193: 160 format (1P,' ',E12.4,' ',E12.4)
195: enddo
196: if (rank .eq. 0) then
197: write(*,*)
198: endif
199: endif
201: ! ** Free work space
202: call EPSDestroy(eps,ierr)
203: call MatDestroy(A,ierr)
204: call VecDestroy(xr,ierr)
205: call VecDestroy(xi,ierr)
207: call SlepcFinalize(ierr)
208: end