MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
DavidsonDiagonalizer_module.f90
Go to the documentation of this file.
1! Copyright 2019
2!
3! For a comprehensive list of the developers that contributed to these codes
4! see the UK-AMOR website.
5!
6! This file is part of UKRmol-in (UKRmol+ suite).
7!
8! UKRmol-in is free software: you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation, either version 3 of the License, or
11! (at your option) any later version.
12!
13! UKRmol-in is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17!
18! You should have received a copy of the GNU General Public License
19! along with UKRmol-in (in source/COPYING). Alternatively, you can also visit
20! <https://www.gnu.org/licenses/>.
21
22!> \brief Diagonalizer type using Davidson backend
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> This type is always available. Uses the SCATCI routine \c mkdvm.
27!>
28!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
29!>
30module davidsondiagonalizer_module
31
32 use precisn, only: wp
33 use const_gbl, only: stdout
34 use baseintegral_module, only: baseintegral
35 use basematrix_module, only: basematrix
36 use writermatrix_module, only: writermatrix
37 use diagonalizer_module, only: basediagonalizer
38 use diagonalizerresult_module, only: diagonalizerresult
39 use options_module, only: options
40 use mpi_gbl, only: master, myrank, mpi_mod_bcast, mpi_xermsg
42
43 implicit none
44
45 type, extends(basediagonalizer) :: davidsondiagonalizer
46 contains
47 procedure, public :: diagonalize => diagonalize_davidson
48 procedure, private :: diagonalize_writermatrix
49 end type davidsondiagonalizer
50
51contains
52
53 subroutine diagonalize_davidson(this, matrix_elements, num_eigenpair, dresult, all_procs, option, integrals)
54 class(davidsondiagonalizer) :: this
55 class(diagonalizerresult) :: dresult
56 class(basematrix), intent(in) :: matrix_elements
57 class(baseintegral), intent(in) :: integrals
58 type(options), intent(in) :: option
59 integer, intent(in) :: num_eigenpair
60 logical, intent(in) :: all_procs
61
62 integer :: max_iterations
63 real(wp) :: max_tolerance
64
65 !Only the first rank can do it
66 if (myrank /= master) return
67
68 max_iterations = option % max_iterations
69 max_tolerance = option % max_tolerance
70
71 if (max_iterations < 0) then
72 max_iterations = max(num_eigenpair * 40, 500)
73 end if
74
75 print *, 'max_iterations', max_iterations
76 print *, 'max_tolerance', max_tolerance
77
78 select type(matrix_elements)
79 type is (writermatrix)
80 call this % diagonalize_writermatrix(matrix_elements, num_eigenpair, dresult, all_procs, &
81 max_iterations, max_tolerance, option, integrals)
82 class is (basematrix)
83 call mpi_xermsg('DavidsonDiagonalizer_module', 'diagonalize_davidson', &
84 'Only WriterMatrix format can use Davidson', 1, 1)
85 end select
86
87 end subroutine diagonalize_davidson
88
89
90 subroutine diagonalize_writermatrix (this, matrix_elements, num_eigenpair, dresult, all_procs, max_iterations, max_tolerance, &
91 option, integrals)
92 class(davidsondiagonalizer) :: this
93 class(diagonalizerresult) :: dresult
94 type(writermatrix), intent(in) :: matrix_elements
95 class(baseintegral), intent(in) :: integrals
96 type(options), intent(in) :: option
97 integer, intent(in) :: num_eigenpair
98 logical, intent(in) :: all_procs
99 integer, intent(in) :: max_iterations
100 real(wp), intent(in) :: max_tolerance
101
102 !for now these are default but will be included i the input
103
104 real(wp) :: CRITC = 1e-10_wp
105 real(wp) :: CRITR = 1e-8_wp
106 real(wp) :: ORTHO = 1e-7_wp
107
108 integer :: matrix_unit
109 integer :: matrix_size
110 integer :: num_matrix_elements_per_record, num_elems
111
112 real(wp), allocatable :: hamiltonian(:), diagonal(:), eig(:)
113 real(wp), pointer :: vec_ptr(:)
114
115 integer :: mpi_num_eig, mpi_num_vecs, ido, ierr
116
117 real(wp), allocatable :: eigen(:)
118 real(wp), allocatable :: vecs(:,:)
119
120 character(len=4), dimension(30) :: NAMP
121 integer, dimension(10) :: NHD
122 integer, dimension(20) :: KEYCSF, NHE
123 real(kind=wp), dimension(41) :: dtnuc
124
125 write (stdout, "('Diagonalization done with Davidson')")
126 write (stdout, "('Parameters:')")
127 write (stdout, "('N: ',i8)") matrix_elements % get_matrix_size()
128 write (stdout, "('Requested # of eigenpairs',i8)") num_eigenpair
129
130 allocate(eigen(num_eigenpair), vecs(matrix_elements % get_matrix_size(), num_eigenpair))
131
132 eigen = 0
133 vecs = 0
134 matrix_size = 0
135
136 if (myrank == master) then
137 !Get our matrix unit
138 matrix_unit = matrix_elements % get_matrix_unit()
139
140 !Lets follow scatcis example even though we can get these information from the matrix class itself
141 rewind matrix_unit
142 read (matrix_unit) matrix_size, num_matrix_elements_per_record, nhd, namp, nhe, dtnuc
143
144 num_elems = matrix_elements % get_size()
145 print *, 'matrix_size', matrix_size
146 print *, 'num_elems', num_elems
147
148 ! allocate(hamiltonian(matrix_size*num_eigenpair),diagonal(matrix_size),eig(matrix_size))
149 allocate(diagonal(matrix_size))
150 call mkdvm(vecs, eigen, diagonal, matrix_size, num_eigenpair, stdout, matrix_unit, &
151 num_matrix_elements_per_record, num_elems, max_tolerance, critc, critr, ortho, &
152 max_iterations, ierr)
153
154 if (ierr /= 0) then
155 call mpi_xermsg('DavidsonDiagonalizer_module', 'DavidsonDiagonalizer', &
156 'Davidson diagonalization failed', ierr, 1)
157 end if
158 end if
159
160 if (all_procs) then
161 call mpi_mod_bcast(matrix_size, master)
162 mpi_num_vecs = matrix_size * num_eigenpair
163 mpi_num_eig = num_eigenpair
164 call mpi_mod_bcast(eigen,master)
165 do ido = 1, num_eigenpair
166 call mpi_mod_bcast(vecs(:,ido), master)
167 end do
168 end if
169
170 if (iand(option % vector_storage_method, pass_to_cdenprop) /= 0) then
171 call dresult % export_header(option, integrals)
172 end if
173 call dresult % write_header(option, integrals)
174 call dresult % handle_eigenvalues(eigen, matrix_elements % diagonal, num_eigenpair, matrix_size)
175
176 do ido = 1, num_eigenpair
177 call dresult % handle_eigenvector(vecs(:,ido), matrix_size)
178 end do
179
180 ! store eigensolutions for the post-processing stage
181 if (iand(dresult % vector_storage, pass_to_cdenprop) /= 0) then
182 !Eigenvectors:
183 if (allocated(dresult % ci_vec % CV)) deallocate(dresult % ci_vec % CV)
184 call move_alloc(vecs, dresult % ci_vec % CV)
185 dresult % ci_vec % CV_is_scalapack = .false.
186
187 !Eigenvalues and phases:
188 call dresult % export_eigenvalues(eigen(1:num_eigenpair), matrix_elements % diagonal, &
189 num_eigenpair, matrix_size, dresult % ci_vec % ei, dresult % ci_vec % iphz)
190
191 write (stdout, '(/," Eigensolutions exported into the CIVect structure.",/)')
192 else
193 deallocate (vecs)
194 end if
195
196 deallocate(eigen)
197
198 end subroutine diagonalize_writermatrix
199
200end module davidsondiagonalizer_module
MPI SCATCI Constants module.
integer, parameter diagonalize
Diagonalize.
integer, parameter pass_to_cdenprop
Eigenpairs will be saved in memory and passed to CDENPROP and outer interface.