Actual source code: ex1f.F
1: !
2: ! Program usage: mpirun -np n ex1f [-help] [-n <n>] [all SLEPc options]
3: !
4: ! Description: Simple example that solves an eigensystem with the EPS object.
5: ! The standard symmetric eigenvalue problem to be solved corresponds to the
6: ! Laplacian operator in 1 dimension.
7: !
8: ! The command line options are:
9: ! -n <n>, where <n> = number of grid points = matrix size
10: !
11: !/*T
12: ! Concepts: SLEPc - Basic functionality
13: ! Routines: SlepcInitialize(); SlepcFinalize();
14: ! Routines: EPSCreate(); EPSSetFromOptions();
15: ! Routines: EPSSolve(); EPSDestroy();
16: !T*/
17: !
18: ! ----------------------------------------------------------------------
19: !
20: program main
21: implicit none
23: #include "finclude/petsc.h"
24: #include "finclude/petscvec.h"
25: #include "finclude/petscmat.h"
26: #include finclude/slepc.h
27: #include finclude/slepceps.h
29: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: ! Declarations
31: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32: !
33: ! Variables:
34: ! A operator matrix
35: ! eps eigenproblem solver context
37: Mat A
38: EPS eps
39: EPSType type
40: PetscReal tol, error
41: PetscScalar kr, ki
42: integer rank, n, nev, ierr, maxit, i, its, nconv
43: integer col(3), Istart, Iend
44: PetscTruth flg
45: PetscScalar value(3)
47: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: ! Beginning of program
49: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
52: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
53: n = 30
54: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
56: if (rank .eq. 0) then
57: write(*,100) n
58: endif
59: 100 format ('1-D Laplacian Eigenproblem, n =',i6)
61: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: ! Compute the operator matrix that defines the eigensystem, Ax=kx
63: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: call MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,A,
66: & ierr)
67: call MatSetFromOptions(A,ierr)
69: call MatGetOwnershipRange(A,Istart,Iend,ierr)
70: if (Istart .eq. 0) then
71: i = 0
72: col(1) = 0
73: col(2) = 1
74: value(1) = 2.0
75: value(2) = -1.0
76: call MatSetValues(A,1,i,2,col,value,INSERT_VALUES,ierr)
77: Istart = Istart+1
78: endif
79: if (Iend .eq. n) then
80: i = n-1
81: col(1) = n-2
82: col(2) = n-1
83: value(1) = -1.0
84: value(2) = 2.0
85: call MatSetValues(A,1,i,2,col,value,INSERT_VALUES,ierr)
86: Iend = Iend-1
87: endif
88: value(1) = -1.0
89: value(2) = 2.0
90: value(3) = -1.0
91: do i=Istart,Iend-1
92: col(1) = i-1
93: col(2) = i
94: col(3) = i+1
95: call MatSetValues(A,1,i,3,col,value,INSERT_VALUES,ierr)
96: enddo
98: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
99: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
101: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: ! Create the eigensolver and display info
103: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: ! ** Create eigensolver context
106: call EPSCreate(PETSC_COMM_WORLD,eps,ierr)
108: ! ** Set operators. In this case, it is a standard eigenvalue problem
109: call EPSSetOperators(eps,A,PETSC_NULL_OBJECT,ierr)
111: ! ** Set solver parameters at runtime
112: call EPSSetFromOptions(eps,ierr)
114: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: ! Solve the eigensystem
116: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: call EPSSolve(eps,ierr)
119: call EPSGetIterationNumber(eps,its,ierr)
120: if (rank .eq. 0) then
121: write(*,*)
122: write(*,140) its
123: endif
124: 140 format (' Number of iterations of the method: ',I4)
125:
126: ! ** Optional: Get some information from the solver and display it
127: call EPSGetType(eps,type,ierr)
128: if (rank .eq. 0) then
129: write(*,110) type
130: endif
131: 110 format (' Solution method: ',A)
132: call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,ierr)
133: if (rank .eq. 0) then
134: write(*,120) nev
135: endif
136: 120 format (' Number of requested eigenvalues:',I2)
137: call EPSGetTolerances(eps,tol,maxit,ierr)
138: if (rank .eq. 0) then
139: write(*,130) tol, maxit
140: endif
141: 130 format (' Stopping condition: tol=',1PE10.4,', maxit=',I6)
143: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: ! Display solution and clean up
145: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: ! ** Get number of converged eigenpairs
148: call EPSGetConverged(eps,nconv,ierr)
149: if (rank .eq. 0) then
150: write(*,150) nconv
151: endif
152: 150 format (' Number of converged approximate eigenpairs:',I2)
154: ! ** Display eigenvalues and relative errors
155: if (nconv.gt.0 .and. rank.eq.0) then
156: write(*,*)
157: write(*,*) ' k ||Ax-kx||/||kx||'
158: write(*,*) ' ----------------- ------------------'
159: do i=0,nconv-1
160: ! ** Get converged eigenpairs: i-th eigenvalue is stored in kr
161: ! ** (real part) and ki (imaginary part)
162: call EPSGetEigenpair(eps,i,kr,ki,PETSC_NULL,PETSC_NULL,ierr)
164: ! ** Compute the relative error associated to each eigenpair
165: call EPSComputeRelativeError(eps,i,error,ierr)
167: if (ki.ne.0.D0) then
168: write(*,180) kr, ki, error
169: else
170: write(*,190) kr, error
171: endif
172: enddo
173: write(*,*)
174: endif
175: 180 format (1P,E11.4,E11.4,' j ',E12.4)
176: 190 format (1P,' ',E12.4,' ',E12.4)
178: ! ** Free work space
179: call EPSDestroy(eps,ierr)
180: call MatDestroy(A,ierr)
182: call SlepcFinalize(ierr)
183: end