MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
Target_RMat_CI_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 Target Rmat CI module
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> Handles the storage of the target coefficients for the contracted Hamiltonian.
27!>
28!> \note 30/01/2017 - Ahmed Al-Refaie: Initial Revision
29!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
30!>
31module target_rmat_ci_module
32
33 use precisn, only: wp
34 use const_gbl, only: stdout
35 use consts_mpi_ci, only: symtype_d2h
36 use scatci_routines, only: cirmat
37 use options_module, only: options, phase
38 use memorymanager_module, only: master_memory
39 use diagonalizerresult_module, only: diagonalizerresult
40 !use BaseCorePotential_module
41
42 implicit none
43
44 public target_rmat_ci, read_ci_mat
45
46 private
47
48 !> \brief This class handles the storage of target CI coefficients.
49 !> \authors A Al-Refaie
50 !> \date 2017
51 !>
52 type, extends(diagonalizerresult) :: target_rmat_ci
53 integer :: set !< Which ci set this belongs to
54 integer :: nstat !< Number of eigenvalues
55 integer :: num_csfs !< Size of the eigenvectors
56 real(wp) :: core_energy !< The associated core energy (for printing)
57 integer :: current_eigenvector = 0
58
59 real(wp), allocatable :: eigenvectors(:,:) !< The ci coefficients
60 real(wp), allocatable :: eigenvalues(:) !< The associated eigenvalues
61 type(phase), pointer :: phase !< Not used, a vestigal type
62
63 !class(BaseCorePotential), pointer :: effective_core_potential => null()
64 contains
65 procedure, public :: initialize => initialize_cirmat
66 procedure, public :: print => print_cirmat
67 procedure, public :: get_coefficient => get_coefficient_ci
68 procedure, public :: print_vecs
69 procedure, public :: export_eigenvalues
70 procedure, public :: handle_eigenvalues => store_eigenvalues
71 procedure, public :: handle_eigenvector => store_eigenvector
72 procedure, public :: modify_diagonal
73 procedure, public :: modify_l2_diagonal
74 end type target_rmat_ci
75
76contains
77
78 !> \brief Sets up and allocates eigenvalues and eigenvectors
79 !> \authors A Al-Refaie
80 !> \date 2017
81 !>
82 !> \param[inout] this Target CI object to initialize.
83 !> \param[in] set Which CI set this belongs to
84 !> \param[in] nstat Number of eigenpairs/target states
85 !> \param[in] num_csfs Size of the eigenvectors/number of ci coefficients per target state
86 !> \param[in] phase_ Not used
87 !> \param[in] core_energy Core energy used for eigenvalue printing
88 !>
89 subroutine initialize_cirmat (this, set, nstat, num_csfs, phase_, core_energy)
90 class(target_rmat_ci) :: this
91 integer, intent(in) :: set, nstat, num_csfs
92 real(wp), intent(in) :: core_energy
93 type(phase), intent(in), target :: phase_
94 !class(BaseCorePotential), pointer :: ecp
95
96 integer :: ifail
97
98 !Set relevant values
99 this % set = set
100 this % nstat = nstat
101 this % num_csfs = num_csfs
102 this % core_energy = core_energy
103 this % phase => phase_
104
105 !if (associated(ecp)) then
106 ! this % effective_core_potential => ecp
107 !end if
108
109 !allocate the necessary arrays and track the memory usage
110 allocate(this % eigenvectors(this % num_csfs, this % nstat), this % eigenvalues(this % nstat), stat = ifail)
111 call master_memory % track_memory(storage_size(this % eigenvectors)/8, &
112 size(this % eigenvectors), ifail, 'TARGET_RMAT_CI::EIGENVECTORS')
113 call master_memory % track_memory(storage_size(this % eigenvalues)/8, &
114 size(this % eigenvalues), ifail, 'TARGET_RMAT_CI::EIGENVALUES')
115
116 !Zero the arrays
117 this % eigenvectors(:,:) = 0.0_wp
118 this % eigenvalues(:) = 0.0_wp
119
120 end subroutine initialize_cirmat
121
122
123 !> \brief Prints the eigen-energies of the target states stored
124 !> \authors A Al-Refaie
125 !> \date 2017
126 !>
127 subroutine print_cirmat (this)
128 class(target_rmat_ci) :: this
129 integer :: I, J
130
131 write (stdout, '(/," SET",I4,4x,A)') &
132 this % set
133 write (stdout, '(/," NOCSF=",I5,4x,"NSTAT=",I5,4x,"EN =",F20.10)') &
134 this % num_csfs, this % nstat, this % core_energy
135 write (stdout, '(/," EIGEN-ENERGIES",/,(16X,5F20.10))') &
136 (this % eigenvalues(i) + this % core_energy, i = 1, this % nstat)
137 write (stdout, '( " STATE",I3,":",I3," transformation vectors selected")') &
138 this % set, this % nstat
139
140 end subroutine print_cirmat
141
142
143 !> \brief Print vectors to standard output
144 !> \authors A Al-Refaie
145 !> \date 2017
146 !>
147 subroutine print_vecs (this)
148 class(target_rmat_ci) :: this
149 integer :: ido, jdo
150
151 do ido = 1, this % nstat
152 write (stdout, *) ' Target state ', this % set, ido
153 do jdo = 1, this % num_csfs
154 write (stdout, *) this % eigenvectors(jdo, ido)
155 end do
156 end do
157
158 end subroutine print_vecs
159
160
161 !> \brief Get the specific coefficient for a target state and configuration
162 !> \authors A Al-Refaie
163 !> \date 2017
164 !>
165 !> \param[in] this Target CI object to query.
166 !> \param[in] target_state Which target state to get from.
167 !> \param[in] configuration_idx Which particular coefficient to get.
168 !>
169 !> \result The resultant coefficient.
170 !>
171 real(wp) function get_coefficient_ci (this, target_state, configuration_idx)
172 class(target_rmat_ci) :: this
173 integer, intent(in) :: target_state, configuration_idx
174
175 get_coefficient_ci = this % eigenvectors(configuration_idx, target_state)
176
177 end function get_coefficient_ci
178
179
180 !> \brief Reads the target coefficients from file
181 !> \authors A Al-Refaie
182 !> \date 2017
183 !>
184 !> \param[in] option An initialized Options object.
185 !> \param[in] cirmats An allocated and initialized array of target CI objects.
186 !>
187 subroutine read_ci_mat (option, cirmats)
188 class(options), intent(in) :: option
189 class(target_rmat_ci), intent(inout) :: cirmats(:)
190
191 real(wp), allocatable :: eigenvectors(:)
192 real(wp), allocatable :: eigenvalues(:)
193
194 integer :: ido, jdo, kdo, ii, jj, ci_size, nstat, ifail
195
196 !Allocate the temporary arrays for eigenvalues and eienvectors
197 allocate(eigenvalues(option % total_num_tgt_states), eigenvectors(option % total_num_ci_target), stat = ifail)
198 call master_memory % track_memory(storage_size(eigenvectors)/8, size(eigenvectors), ifail, &
199 'TARGET_RMAT_CI::READ::EIGENVECTORS')
200 call master_memory % track_memory(storage_size(eigenvalues)/8, size(eigenvalues), ifail, &
201 'TARGET_RMAT_CI::READ::EIGENVALUES')
202
203 !Use SCATCI to read them
204 call cirmat(option % num_target_sym, &
205 option % ci_target_switch, &
206 option % ci_set_number, &
207 eigenvectors, &
208 eigenvalues, &
209 option % ci_sequence_number, &
210 option % num_ci_target_sym, &
211 stdout, &
212 option % num_target_state_sym, &
213 option % phase_index, &
214 option % num_continuum, &
215 option % sym_group_flag)
216
217 ii = 0
218 jj = 0
219
220 !Seperate out all the eigenvectors and values into distinct target symmetries and states.
221 do ido = 1, option % num_target_sym
222 ci_size = cirmats(ido) % num_csfs
223 nstat = cirmats(ido) % nstat
224 do jdo = 1, nstat
225 jj = jj + 1
226 cirmats(ido) % eigenvalues(jdo) = eigenvalues(jj)
227 cirmats(ido) % eigenvectors(1:ci_size,jdo) = eigenvectors(ii+1:ii+ci_size)
228 ii = ii + ci_size
229 end do
230 end do
231
232 !Stop tracking this memory
233 call master_memory % free_memory(storage_size(eigenvectors)/8, size(eigenvectors))
234 call master_memory % free_memory(storage_size(eigenvalues)/8, size(eigenvalues))
235
236 !Cleanup
237 deallocate(eigenvalues)
238 deallocate(eigenvectors)
239
240 end subroutine read_ci_mat
241
242
243 !> \brief Save eigenvalues from diagonalizer
244 !> \authors A Al-Refaie
245 !> \date 2017
246 !>
247 !> Copy eigenvalues obtained by the diagonalizer to internal storage of Target_RMat_CI.
248 !>
249 subroutine store_eigenvalues (this, eigenvalues, diagonals, num_eigenpairs, vec_dimen)
250 class(target_rmat_ci) :: this
251 integer, intent(in) :: num_eigenpairs, vec_dimen
252 real(wp), intent(in) :: eigenvalues(num_eigenpairs), diagonals(vec_dimen)
253
254 this % eigenvalues = eigenvalues
255
256 !if (associated(this % effective_core_potential)) then
257 ! call this % effective_core_potential % modify_target_energies(this % set, this % eigenvalues)
258 !end if
259
260 end subroutine store_eigenvalues
261
262
263 !> \brief Save eigenvector from diagonalizer
264 !> \authors A Al-Refaie
265 !> \date 2017
266 !>
267 !> Copy eigenvector obtained by the diagonalizer to internal storage of Target_RMat_CI.
268 !>
269 subroutine store_eigenvector (this, eigenvector, vec_dimen)
270 class(target_rmat_ci) :: this
271 integer, intent(in) :: vec_dimen
272 real(wp), intent(in) :: eigenvector(vec_dimen)
273
274 this % current_eigenvector = this % current_eigenvector + 1
275 this % eigenvectors(:, this % current_eigenvector) = eigenvector(:)
276
277 end subroutine store_eigenvector
278
279
280 !> \brief To be overriden by derived types
281 !> \authors J Benda
282 !> \date 2018
283 !>
284 !> This subroutine body was added so that Target_RMat_CI does not need to be `abstract`.
285 !> Otherwise there were problems with compilation using non-Intel compilers.
286 !>
287 subroutine export_eigenvalues (this, eigenvalues, diagonals, num_eigenpairs, vec_dimen, ei, iphz)
288 use precisn, only: wp
289 class(target_rmat_ci) :: this
290 integer, intent(in) :: num_eigenpairs, vec_dimen
291 real(wp), intent(in) :: eigenvalues(num_eigenpairs), diagonals(vec_dimen)
292 real(wp), allocatable :: ei(:)
293 integer, allocatable :: iphz(:)
294 end subroutine export_eigenvalues
295
296
297 !> \brief Dummy subroutine
298 !> \authors A Al-Refaie
299 !> \date 2017
300 !>
301 !> Has no effect.
302 !>
303 subroutine modify_diagonal (this, value)
304 class(target_rmat_ci), intent(in) :: this
305 real(wp), intent(inout) :: value
306
307 !if (associated(this % effective_core_potential)) then
308 ! call this % effective_core_potential % modify_ecp_diagonal(value)
309 !end if
310
311 end subroutine modify_diagonal
312
313
314 !> \brief Dummy subroutine
315 !> \authors A Al-Refaie
316 !> \date 2017
317 !>
318 !> Has no effect.
319 !>
320 subroutine modify_l2_diagonal (this, value)
321 class(target_rmat_ci), intent(in) :: this
322 real(wp), intent(inout) :: value
323
324 !if (associated(this % effective_core_potential)) then
325 ! call this % effective_core_potential % modify_ecp_diagonal_L2(value)
326 !end if
327
328 end subroutine modify_l2_diagonal
329
330end module target_rmat_ci_module
MPI SCATCI Constants module.
integer, parameter symtype_d2h
This describes D_2_h symmetries.