MPI-SCATCI  2.0
An MPI version of SCATCI
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 
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
40  !use BaseCorePotential_module
41 
42  implicit none
43 
45 
46  private
47 
53  integer :: set
54  integer :: nstat
55  integer :: num_csfs
56  real(wp) :: core_energy
57  integer :: current_eigenvector = 0
58 
59  real(wp), allocatable :: eigenvectors(:,:)
60  real(wp), allocatable :: eigenvalues(:)
61  type(phase), pointer :: phase
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 
76 contains
77 
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(kind(this % eigenvectors), &
112  size(this % eigenvectors), ifail, 'TARGET_RMAT_CI::EIGENVECTORS')
113  call master_memory % track_memory(kind(this % eigenvalues), &
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 
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 
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 
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 
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(kind(eigenvectors), size(eigenvectors), ifail, 'TARGET_RMAT_CI::READ::EIGENVECTORS')
199  call master_memory % track_memory(kind(eigenvalues), size(eigenvalues), ifail, 'TARGET_RMAT_CI::READ::EIGENVALUES')
200 
201  !Use SCATCI to read them
202  call cirmat(option % num_target_sym, &
203  option % ci_target_switch, &
204  option % ci_set_number, &
205  eigenvectors, &
206  eigenvalues, &
207  option % ci_sequence_number, &
208  option % num_ci_target_sym, &
209  stdout, &
210  option % num_target_state_sym, &
211  option % phase_index, &
212  option % num_continuum, &
213  option % sym_group_flag)
214 
215  ii = 0
216  jj = 0
217 
218  !Seperate out all the eigenvectors and values into distinct target symmetries and states.
219  do ido = 1, option % num_target_sym
220  ci_size = cirmats(ido) % num_csfs
221  nstat = cirmats(ido) % nstat
222  do jdo = 1, nstat
223  jj = jj + 1
224  cirmats(ido) % eigenvalues(jdo) = eigenvalues(jj)
225  cirmats(ido) % eigenvectors(1:ci_size,jdo) = eigenvectors(ii+1:ii+ci_size)
226  ii = ii + ci_size
227  end do
228  end do
229 
230  !Stop tracking this memory
231  call master_memory % free_memory(kind(eigenvectors), size(eigenvectors))
232  call master_memory % free_memory(kind(eigenvalues), size(eigenvalues))
233 
234  !Cleanup
235  deallocate(eigenvalues)
236  deallocate(eigenvectors)
237 
238  end subroutine read_ci_mat
239 
240 
247  subroutine store_eigenvalues (this, eigenvalues, diagonals, num_eigenpairs, vec_dimen)
248  class(target_rmat_ci) :: this
249  integer, intent(in) :: num_eigenpairs, vec_dimen
250  real(wp), intent(in) :: eigenvalues(num_eigenpairs), diagonals(vec_dimen)
251 
252  this % eigenvalues = eigenvalues
253 
254  !if (associated(this % effective_core_potential)) then
255  ! call this % effective_core_potential % modify_target_energies(this % set, this % eigenvalues)
256  !end if
257 
258  end subroutine store_eigenvalues
259 
260 
267  subroutine store_eigenvector (this, eigenvector, vec_dimen)
268  class(target_rmat_ci) :: this
269  integer, intent(in) :: vec_dimen
270  real(wp), intent(in) :: eigenvector(vec_dimen)
271 
272  this % current_eigenvector = this % current_eigenvector + 1
273  this % eigenvectors(:, this % current_eigenvector) = eigenvector(:)
274 
275  end subroutine store_eigenvector
276 
277 
285  subroutine export_eigenvalues (this, eigenvalues, diagonals, num_eigenpairs, vec_dimen, ei, iphz)
286  use precisn, only: wp
287  class(target_rmat_ci) :: this
288  integer, intent(in) :: num_eigenpairs, vec_dimen
289  real(wp), intent(in) :: eigenvalues(num_eigenpairs), diagonals(vec_dimen)
290  real(wp), allocatable :: ei(:)
291  integer, allocatable :: iphz(:)
292  end subroutine export_eigenvalues
293 
294 
301  subroutine modify_diagonal (this, value)
302  class(target_rmat_ci), intent(in) :: this
303  real(wp), intent(inout) :: value
304 
305  !if (associated(this % effective_core_potential)) then
306  ! call this % effective_core_potential % modify_ecp_diagonal(value)
307  !end if
308 
309  end subroutine modify_diagonal
310 
311 
318  subroutine modify_l2_diagonal (this, value)
319  class(target_rmat_ci), intent(in) :: this
320  real(wp), intent(inout) :: value
321 
322  !if (associated(this % effective_core_potential)) then
323  ! call this % effective_core_potential % modify_ecp_diagonal_L2(value)
324  !end if
325 
326  end subroutine modify_l2_diagonal
327 
328 end module target_rmat_ci_module
diagonalizerresult_module::diagonalizerresult
Output from diagonalization.
Definition: DiagonalizerResult_module.f90:48
memorymanager_module
Memory manager module.
Definition: MemoryManager_module.f90:31
target_rmat_ci_module::get_coefficient_ci
real(wp) function get_coefficient_ci(this, target_state, configuration_idx)
Get the specific coefficient for a target state and configuration.
Definition: Target_RMat_CI_module.f90:172
diagonalizerresult_module::handle_eigenvector
Main build routine of the hamiltonian.
Definition: DiagonalizerResult_module.f90:87
consts_mpi_ci::symtype_d2h
integer, parameter symtype_d2h
This describes D_2_h symmetries.
Definition: consts_mpi_ci.f90:57
target_rmat_ci_module::target_rmat_ci
This class handles the storage of target CI coefficients.
Definition: Target_RMat_CI_module.f90:52
options_module::options
Calculation setup for a single Hamiltonian diagonalization.
Definition: Options_module.f90:74
target_rmat_ci_module::export_eigenvalues
subroutine export_eigenvalues(this, eigenvalues, diagonals, num_eigenpairs, vec_dimen, ei, iphz)
To be overriden by derived types.
Definition: Target_RMat_CI_module.f90:286
target_rmat_ci_module::store_eigenvector
subroutine store_eigenvector(this, eigenvector, vec_dimen)
Save eigenvector from diagonalizer.
Definition: Target_RMat_CI_module.f90:268
target_rmat_ci_module::read_ci_mat
subroutine, public read_ci_mat(option, cirmats)
Reads the target coefficients from file.
Definition: Target_RMat_CI_module.f90:188
consts_mpi_ci
MPI SCATCI Constants module.
Definition: consts_mpi_ci.f90:31
target_rmat_ci_module::store_eigenvalues
subroutine store_eigenvalues(this, eigenvalues, diagonals, num_eigenpairs, vec_dimen)
Save eigenvalues from diagonalizer.
Definition: Target_RMat_CI_module.f90:248
target_rmat_ci_module::modify_diagonal
subroutine modify_diagonal(this, value)
Dummy subroutine.
Definition: Target_RMat_CI_module.f90:302
diagonalizerresult_module::export_eigenvalues
build routine of the hamiltonian
Definition: DiagonalizerResult_module.f90:103
options_module::phase
?
Definition: Options_module.f90:61
diagonalizerresult_module
Diagonalizer result data object.
Definition: DiagonalizerResult_module.f90:28
target_rmat_ci_module::print_cirmat
subroutine print_cirmat(this)
Prints the eigen-energies of the target states stored.
Definition: Target_RMat_CI_module.f90:128
target_rmat_ci_module::modify_l2_diagonal
subroutine modify_l2_diagonal(this, value)
Dummy subroutine.
Definition: Target_RMat_CI_module.f90:319
target_rmat_ci_module::initialize_cirmat
subroutine initialize_cirmat(this, set, nstat, num_csfs, phase_, core_energy)
Sets up and allocates eigenvalues and eigenvectors.
Definition: Target_RMat_CI_module.f90:90
options_module
Options module.
Definition: Options_module.f90:41
target_rmat_ci_module::print_vecs
subroutine print_vecs(this)
Print vectors to standard output.
Definition: Target_RMat_CI_module.f90:148
diagonalizerresult_module::handle_eigenvalues
Main build routine of the hamiltonian.
Definition: DiagonalizerResult_module.f90:71
target_rmat_ci_module
Target Rmat CI module.
Definition: Target_RMat_CI_module.f90:31
memorymanager_module::master_memory
type(memorymanager), public master_memory
This is the global memory manager.
Definition: MemoryManager_module.f90:69