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: !/*T
13: ! Concepts: SLEPc - Basic functionality
14: ! Routines: SlepcInitialize(); SlepcFinalize();
15: ! Routines: EPSCreate(); EPSSetFromOptions();
16: ! Routines: EPSSolve(); EPSDestroy();
17: !T*/
18: !
19: ! ----------------------------------------------------------------------
20: !
21: program main
22: implicit none
24: #include "finclude/petsc.h"
25: #include "finclude/petscvec.h"
26: #include "finclude/petscmat.h"
27: #include finclude/slepc.h
28: #include finclude/slepceps.h
30: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: ! Declarations
32: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: !
34: ! Variables:
35: ! A operator matrix
36: ! eps eigenproblem solver context
38: Mat A
39: EPS eps
40: EPSType type
41: PetscReal tol, error
42: PetscScalar kr, ki
43: integer size, rank, N, m, nev, ierr, maxit, i, its, nconv
44: PetscTruth flg
46: ! This is the routine to use for matrix-free approach
47: !
48: external MatIsing_Mult
50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: ! Beginning of program
52: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
55: #if defined(PETSC_USE_COMPLEX)
56: write(*,*) 'This example requires real numbers.'
57: goto 999
58: #endif
59: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
60: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
61: if (size .ne. 1) then
62: if (rank .eq. 0) then
63: write(*,*) 'This is a uniprocessor example only!'
64: endif
65: SETERRQ(1,' ',ierr)
66: endif
67: m = 30
68: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-m',m,flg,ierr)
69: N = 2*m
71: if (rank .eq. 0) then
72: write(*,*)
73: write(*,'(A,I6,A)') 'Ising Model Eigenproblem, m=',m,', (N=2*m)'
74: write(*,*)
75: endif
77: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: ! Register the matrix-vector subroutine for the operator that defines
79: ! the eigensystem, Ax=kx
80: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: call MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,PETSC_NULL_OBJECT,A,
83: & ierr)
84: call MatShellSetOperation(A,MATOP_MULT,MatIsing_Mult,ierr)
86: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: ! Create the eigensolver and display info
88: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: ! ** Create eigensolver context
91: call EPSCreate(PETSC_COMM_WORLD,eps,ierr)
93: ! ** Set operators. In this case, it is a standard eigenvalue problem
94: call EPSSetOperators(eps,A,PETSC_NULL_OBJECT,ierr)
95: call EPSSetProblemType(eps,EPS_NHEP,ierr)
97: ! ** Set solver parameters at runtime
98: call EPSSetFromOptions(eps,ierr)
100: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: ! Solve the eigensystem
102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: call EPSSolve(eps,ierr)
105: call EPSGetIterationNumber(eps,its,ierr)
106: if (rank .eq. 0) then
107: write(*,'(A,I4)') ' Number of iterations of the method: ', its
108: endif
110: ! ** Optional: Get some information from the solver and display it
111: call EPSGetType(eps,type,ierr)
112: if (rank .eq. 0) then
113: write(*,'(A,A)') ' Solution method: ', type
114: endif
115: call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,ierr)
116: if (rank .eq. 0) then
117: write(*,'(A,I2)') ' Number of requested eigenvalues:', nev
118: endif
119: call EPSGetTolerances(eps,tol,maxit,ierr)
120: if (rank .eq. 0) then
121: write(*,'(A,1PE10.4,A,I6)') ' Stopping condition: tol=', tol,
122: & ', maxit=', maxit
123: endif
125: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: ! Display solution and clean up
127: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: ! ** Get number of converged eigenpairs
130: call EPSGetConverged(eps,nconv,ierr)
131: if (rank .eq. 0) then
132: write(*,'(A,I2)') ' Number of converged eigenpairs:', nconv
133: endif
135: ! ** Display eigenvalues and relative errors
136: if (nconv.gt.0 .and. rank.eq.0) then
137: write(*,*)
138: write(*,*) ' k ||Ax-kx||/||kx||'
139: write(*,*) ' ----------------- ------------------'
140: do i=0,nconv-1
141: ! ** Get converged eigenpairs: i-th eigenvalue is stored in kr
142: ! ** (real part) and ki (imaginary part)
143: call EPSGetEigenpair(eps,i,kr,ki,PETSC_NULL,PETSC_NULL,ierr)
145: ! ** Compute the relative error associated to each eigenpair
146: call EPSComputeRelativeError(eps,i,error,ierr)
148: if (ki.ne.0.D0) then
149: write(*,'(1P,E11.4,E11.4,A,E12.4)') kr, ki, ' j ', error
150: else
151: write(*,'(1P,A,E12.4,A,E12.4)') ' ', kr, ' ', error
152: endif
153: enddo
154: write(*,*)
155: endif
157: ! ** Free work space
158: call EPSDestroy(eps,ierr)
159: call MatDestroy(A,ierr)
161: 999 continue
162: call SlepcFinalize(ierr)
163: end
165: ! -------------------------------------------------------------------
166: !
167: ! MatIsing_Mult - user provided matrix-vector multiply
168: !
169: ! Input Parameters:
170: ! A - matrix
171: ! x - input vector
172: !
173: ! Output Parameter:
174: ! y - output vector
175: !
176: subroutine MatIsing_Mult(A,x,y,ierr)
177: implicit none
179: #include "finclude/petsc.h"
181: Mat A
182: Vec x,y
183: integer trans,one,ierr,i,N
184: PetscScalar x_array(1),y_array(1)
185: PetscOffset i_x,i_y
187: ! The actual routine for the matrix-vector product
188: external mvmisg
190: call MatGetSize(A,N,PETSC_NULL_INTEGER,ierr)
191: call VecGetArray(x,x_array,i_x,ierr)
192: call VecGetArray(y,y_array,i_y,ierr)
194: trans = 0
195: one = 1
196: call mvmisg(trans,N,one,x_array(i_x+1),N,y_array(i_y+1),N)
198: call VecRestoreArray(x,x_array,i_x,ierr)
199: call VecRestoreArray(y,y_array,i_y,ierr)
201: return
202: end