32 use const_gbl,
only: stdout
34 use mpi_gbl,
only: master
35 use precisn,
only: longint, wp
49 #include <finclude/petsc.h>
50 #include <finclude/slepceps.h>
71 mat,
pointer,
intent(inout) :: mat_out
73 write(stdout,
"('--------Optimized SLEPC Matrix Format chosen------')")
75 mat_out => matrix_elements % get_PETSC_matrix()
77 if (.not.
associated(mat_out))
then
78 stop
"Matrix doesn;t exist for SLEPC, are you sure you created it?"
86 eps,
intent(inout) :: eps
87 petscerrorcode :: ierr
89 write(stdout,
"('--------KRYLOVSCHUR used as Diagonalizer------')")
92 call epssettype(eps,epskrylovschur,ierr)
100 subroutine diagonalize_slepc (this, matrix_elements, num_eigenpair, dresult, all_procs, option, integrals)
103 class(
basematrix),
intent(in) :: matrix_elements
105 type(
options),
intent(in) :: option
106 integer,
intent(in) :: num_eigenpair
107 logical,
intent(in) :: all_procs
109 integer :: mat_dimen, number_solved
110 real(wp),
allocatable :: diagonal_temp(:), eigen(:)
117 epsconvergedreason :: conv_reason
118 petscreal :: maxtol,tol, error
119 petscscalar :: kr, ki, current_shift
120 petscscalar,
pointer :: xx_v(:)
122 petscint :: nev, maxit, its, nconv, matrix_size, ido, n
123 petscerrorcode :: ierr
128 maxit = option % max_iterations
129 maxtol = option % max_tolerance
131 if (maxit < 0) maxit = petsc_default_integer
133 matrix_size = matrix_elements % get_matrix_size()
134 mat_dimen = matrix_size
137 write (stdout,
"('-----------WARNING----------------------- ')")
138 write (stdout,
"('When dealing with symmetries that are degenerate SLEPC may give incorrect values')")
139 write (stdout,
"('YOU HAVE BEEN WARNED')")
140 write (stdout,
"('---------------------------------------- ')")
141 write (stdout,
"('N: ',i8)") matrix_size
142 write (stdout,
"('Requested # of eigenpairs',i8)") num_eigenpair
143 write (stdout,
"('Maximum iterations',i8)") maxit
144 write (stdout,
"('Maximum tolerance',es9.2)") maxtol
147 select type(matrix_elements)
149 call this % UseSLEPCMat(matrix_elements, a)
151 stop
"Matrix format not yet implemented for SLEPC"
156 call matcreatevecs(a, xr, xi, ierr)
159 call vecscattercreatetoall(xr, vsctx, f_xr, ierr)
161 call vecscattercreatetozero(xr, vsctx, f_xr, ierr)
164 call epscreate(int(grid % gcomm, kind(petsc_comm_world)), eps, ierr)
165 call this % SelectEPS(eps)
166 call epssetoperators(eps, a, petsc_null_mat, ierr)
167 call epssetproblemtype(eps, eps_hep, ierr)
200 call epssetdimensions(eps, n, petsc_decide, petsc_decide, ierr)
202 call epssettolerances(eps, maxtol, maxit, ierr)
205 call epssolve(eps, ierr)
211 call epsgetiterationnumber(eps, its, ierr)
212 call epsgettolerances(eps, tol, maxit, ierr)
214 if (grid % grank == master)
then
216 write (stdout,
"(('/ Stopping condition: tol=',1P,E11.4,', maxit=',I4,', stopped at it=',I4))") tol, maxit, its
219 call epsgetconverged(eps, nconv, ierr)
221 if (nconv < num_eigenpair)
then
222 write (stdout, *)
'Not all requested eigenpairs have converged!!!!'
223 write (stdout, *)
'Only ', nconv,
' have converged against ', num_eigenpair,
' requested'
224 stop
"Diagoanlization not completed"
227 allocate(eigen(num_eigenpair))
229 do ido = 0, min(nconv, n) - 1
232 call epsgeteigenpair(eps, ido, kr, ki, xr, xi, ierr)
233 eigen(ido + 1) = petscrealpart(kr)
238 call dresult % export_header(option, integrals)
242 if (grid % grank == master)
then
243 call dresult % write_header(option, integrals)
247 if (all_procs .or. grid % grank == master)
then
248 call dresult % handle_eigenvalues(eigen, matrix_elements % diagonal, num_eigenpair, mat_dimen)
251 do ido = 0, min(nconv, n) - 1
254 call epsgeteigenpair(eps, ido, kr, ki, xr, xi, ierr)
264 call vecscatterbegin(vsctx, xr, f_xr, insert_values, scatter_forward, ierr)
265 call vecscatterend(vsctx, xr, f_xr, insert_values, scatter_forward, ierr)
266 if (all_procs .or. grid % grank == master)
then
267 call vecgetarrayreadf90(f_xr, xx_v, ierr)
268 call dresult % handle_eigenvector(xx_v(1:mat_dimen), mat_dimen)
270 call vecrestorearrayreadf90(f_xr, xx_v, ierr)
275 if (grid % grank == master)
then
276 call dresult % finalize_solutions(option)
279 call epsdestroy(eps, ierr)
280 call vecdestroy(xr, ierr)
281 call vecdestroy(xi, ierr)