Actual source code: ex6f.F
1: !
2: ! Program usage: mpirun -np n ex6f [-help] [-m <m>] [all SLEPc options]
3: !
4: ! Description: This example solves the eigensystem arising in the Ising
5: ! model for ferromagnetic materials. The file mvmisg.f must be linked
6: ! together. Information about the model can be found at the following
7: ! site http://math.nist.gov/MatrixMarket/data/NEP
8: !
9: ! The command line options are:
10: ! -m <m>, where <m> is the number of 2x2 blocks, i.e. matrix size N=2*m
11: !
12: ! ----------------------------------------------------------------------
13: !
14: program main
15: implicit none
17: #include "finclude/petsc.h"
18: #include "finclude/petscvec.h"
19: #include "finclude/petscmat.h"
20: #include finclude/slepc.h
21: #include finclude/slepceps.h
23: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: ! Declarations
25: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26: !
27: ! Variables:
28: ! A operator matrix
29: ! eps eigenproblem solver context
31: Mat A
32: EPS eps
33: EPSType type
34: PetscReal tol, error
35: PetscScalar kr, ki
36: integer size, rank, N, m, nev, ierr, maxit, i, its, nconv
37: PetscTruth flg
39: ! This is the routine to use for matrix-free approach
40: !
41: external MatIsing_Mult
43: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: ! Beginning of program
45: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
48: #if defined(PETSC_USE_COMPLEX)
49: write(*,*) 'This example requires real numbers.'
50: goto 999
51: #endif
52: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
53: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
54: if (size .ne. 1) then
55: if (rank .eq. 0) then
56: write(*,*) 'This is a uniprocessor example only!'
57: endif
58: SETERRQ(1,' ',ierr)
59: endif
60: m = 30
61: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-m',m,flg,ierr)
62: N = 2*m
64: if (rank .eq. 0) then
65: write(*,*)
66: write(*,'(A,I6,A)') 'Ising Model Eigenproblem, m=',m,', (N=2*m)'
67: write(*,*)
68: endif
70: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: ! Register the matrix-vector subroutine for the operator that defines
72: ! the eigensystem, Ax=kx
73: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: call MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,PETSC_NULL_OBJECT,A,
76: & ierr)
77: call MatShellSetOperation(A,MATOP_MULT,MatIsing_Mult,ierr)
79: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: ! Create the eigensolver and display info
81: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: ! ** Create eigensolver context
84: call EPSCreate(PETSC_COMM_WORLD,eps,ierr)
86: ! ** Set operators. In this case, it is a standard eigenvalue problem
87: call EPSSetOperators(eps,A,PETSC_NULL_OBJECT,ierr)
88: call EPSSetProblemType(eps,EPS_NHEP,ierr)
90: ! ** Set solver parameters at runtime
91: call EPSSetFromOptions(eps,ierr)
93: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: ! Solve the eigensystem
95: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: call EPSSolve(eps,ierr)
98: call EPSGetIterationNumber(eps,its,ierr)
99: if (rank .eq. 0) then
100: write(*,'(A,I4)') ' Number of iterations of the method: ', its
101: endif
103: ! ** Optional: Get some information from the solver and display it
104: call EPSGetType(eps,type,ierr)
105: if (rank .eq. 0) then
106: write(*,'(A,A)') ' Solution method: ', type
107: endif
108: call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,ierr)
109: if (rank .eq. 0) then
110: write(*,'(A,I2)') ' Number of requested eigenvalues:', nev
111: endif
112: call EPSGetTolerances(eps,tol,maxit,ierr)
113: if (rank .eq. 0) then
114: write(*,'(A,1PE10.4,A,I6)') ' Stopping condition: tol=', tol,
115: & ', maxit=', maxit
116: endif
118: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119: ! Display solution and clean up
120: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: ! ** Get number of converged eigenpairs
123: call EPSGetConverged(eps,nconv,ierr)
124: if (rank .eq. 0) then
125: write(*,'(A,I2)') ' Number of converged eigenpairs:', nconv
126: endif
128: ! ** Display eigenvalues and relative errors
129: if (nconv.gt.0 .and. rank.eq.0) then
130: write(*,*)
131: write(*,*) ' k ||Ax-kx||/||kx||'
132: write(*,*) ' ----------------- ------------------'
133: do i=0,nconv-1
134: ! ** Get converged eigenpairs: i-th eigenvalue is stored in kr
135: ! ** (real part) and ki (imaginary part)
136: call EPSGetEigenpair(eps,i,kr,ki,PETSC_NULL,PETSC_NULL,ierr)
138: ! ** Compute the relative error associated to each eigenpair
139: call EPSComputeRelativeError(eps,i,error,ierr)
141: if (ki.ne.0.D0) then
142: write(*,'(1P,E11.4,E11.4,A,E12.4)') kr, ki, ' j ', error
143: else
144: write(*,'(1P,A,E12.4,A,E12.4)') ' ', kr, ' ', error
145: endif
146: enddo
147: endif
148: write(*,*)
150: ! ** Free work space
151: call EPSDestroy(eps,ierr)
152: call MatDestroy(A,ierr)
154: 999 continue
155: call SlepcFinalize(ierr)
156: end
158: ! -------------------------------------------------------------------
159: !
160: ! MatIsing_Mult - user provided matrix-vector multiply
161: !
162: ! Input Parameters:
163: ! A - matrix
164: ! x - input vector
165: !
166: ! Output Parameter:
167: ! y - output vector
168: !
169: subroutine MatIsing_Mult(A,x,y,ierr)
170: implicit none
172: #include "finclude/petsc.h"
174: Mat A
175: Vec x,y
176: integer trans,one,ierr,N
177: PetscScalar x_array(1),y_array(1)
178: PetscOffset i_x,i_y
180: ! The actual routine for the matrix-vector product
181: external mvmisg
183: call MatGetSize(A,N,PETSC_NULL_INTEGER,ierr)
184: call VecGetArray(x,x_array,i_x,ierr)
185: call VecGetArray(y,y_array,i_y,ierr)
187: trans = 0
188: one = 1
189: call mvmisg(trans,N,one,x_array(i_x+1),N,y_array(i_y+1),N)
191: call VecRestoreArray(x,x_array,i_x,ierr)
192: call VecRestoreArray(y,y_array,i_y,ierr)
194: return
195: end