MPI-SCATCI  2.0
An MPI version of SCATCI
/home/jakub/Dokumenty/codes/ukrmol-in-git/source/mpi-ci-diag/mpi_scatci.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 
41 program mpi_scatci
42 
44  use const_gbl, only: stdout
45  use global_utils, only: print_ukrmol_header
46  use precisn, only: wp, longint
47  use mpi_gbl, only: mpi_mod_start, mpi_mod_finalize, master, myrank, mpi_mod_print_info
49  use orbital_module, only: orbitaltable
52  use basematrix_module, only: basematrix
60  use timing_module, only: master_timer
65 
66  implicit none
67 
68  type(optionsset) :: scatci_input
69  type(orbitaltable) :: orbitals
70  type(uncontracted_hamiltonian) :: enrgms
71  type(contracted_hamiltonian) :: enrgmx
72  type(csfmanager) :: configuration_manager
73  class(baseintegral), pointer :: integrals => null()
74  class(target_ci_hamiltonian), pointer :: tgt_ci_hamiltonian => null()
75  class(basediagonalizer), pointer :: diagonalizer => null()
76  class(basematrix), pointer :: matrix_elements => null()
77  type(solutionhandler), allocatable :: solutions(:)
78  type(csfobject), allocatable :: csfs(:)
79  type(target_rmat_ci), allocatable :: ci_rmat(:)
80  real(wp), allocatable :: test_eig(:), test_vecs(:,:)
81  integer :: num_mat_elms, i, j
82 
83  logical, parameter :: sequential_diagonalizations = .true.
84  logical, parameter :: master_writes_to_stdout = .true.
85  logical, parameter :: allow_shared_memory = .true.
86 
87  call mpi_mod_start(master_writes_to_stdout, allow_shared_memory)
88  call master_timer % initialize()
89  if (myrank == master) call print_ukrmol_header(stdout)
90  call mpi_mod_print_info(stdout)
91 
92  ! Initialize all libraries used by MPI-SCATCI that need it (e.g. SLEPc)
94 
95  ! Read options (for all symmetries)
96  call scatci_input % read
97  allocate (solutions(size(scatci_input % opts)))
98 
99  ! Set up the MPI/BLACS groups
100  call process_grid % setup(size(scatci_input % opts), sequential_diagonalizations)
101  call process_grid % summarize
102  if (.not. process_grid % sequential) call scatci_input % setup_write_order
103 
104  ! Process all symmetries
105  symmetry_loop: do i = 1, size(scatci_input % opts)
106 
107  if (.not. process_grid % is_my_group_work(i)) cycle
108 
109  call master_timer % start_timer("Wall Time")
110 
111  call orbitals % initialize(scatci_input % opts(i) % tot_num_spin_orbitals, &
112  scatci_input % opts(i) % tot_num_orbitals, &
113  scatci_input % opts(i) % sym_group_flag, &
114  scatci_input % opts(i) % num_syms, &
115  scatci_input % opts(i) % positron_flag)
116  call orbitals % construct(scatci_input % opts(i) % num_orbitals_sym, &
117  scatci_input % opts(i) % num_target_orbitals_sym_dinf_congen, &
118  scatci_input % opts(i) % num_orbitals_sym, &
119  scatci_input % opts(i) % num_electronic_orbitals_sym, &
120  scatci_input % opts(i) % num_target_orbitals_sym_congen)
121 
122  ! Compute the electron number for each reference orbital
123  call orbitals % compute_electron_index(scatci_input % opts(i) % num_electrons, &
124  scatci_input % opts(i) % reference_dtrs)
125 
126  ! Construct our csf manager
127  call configuration_manager % initialize(scatci_input % opts(i), orbitals)
128 
129  ! Construct our CSFS
130  call master_timer % start_timer("Construct CSFs")
131  call configuration_manager % create_csfs(csfs, &
132  scatci_input % opts(i) % orbital_sequence_number, &
133  scatci_input % opts(i) % num_ci_target_sym, &
134  scatci_input % opts(i) % lambda_continuum_orbitals_target)
135  call master_timer % stop_timer("Construct CSFs")
136  call master_memory % print_memory_report
137 
138  ! Read molecular integrals (assuming all symmetries use the same!)
139  if (.not. associated(integrals)) then
140  call dispatchintegral(scatci_input % opts(i) % sym_group_flag, &
141  scatci_input % opts(i) % use_UKRMOL_integrals, &
142  integrals)
143  call integrals % initialize(scatci_input % opts(i), &
144  orbitals % total_num_orbitals, &
145  orbitals % orbital_map)
146  call integrals % load_integrals(scatci_input % opts(i) % integral_unit)
147  call master_memory % print_memory_report
148  end if
149 
150  write (stdout, *) ' Load CI target and construct matrix elements'
151 
152  ! Before we do anything else let us get our ECP if we need them
153  !call DispatchECP(SCATCI_input % opts(i) % ecp_type, &
154  ! SCATCI_input % opts(i) % ecp_filename, &
155  ! SCATCI_input % opts(i) % all_ecp_defined, &
156  ! SCATCI_input % opts(i) % num_target_sym, &
157  ! SCATCI_input % opts(i) % num_target_state_sym, &
158  ! SCATCI_input % opts(i) % target_spatial, &
159  ! SCATCI_input % opts(i) % target_multiplicity, ecp)
160 
161  if (scatci_input % opts(i) % do_ci_contraction()) then
162  ! allocate the CIRmats
163  allocate(ci_rmat(scatci_input % opts(i) % num_target_sym))
164  do j = 1, scatci_input % opts(i) % num_target_sym
165  call ci_rmat(j) % initialize(j, &
166  scatci_input % opts(i) % num_target_state_sym(j), &
167  scatci_input % opts(i) % num_ci_target_sym(j), &
168  scatci_input % opts(i) % ci_phase(j), &
169  integrals % get_core_energy())
170  end do
171 
172  if (scatci_input % opts(i) % ci_target_switch > 0) then
173  call read_ci_mat(scatci_input % opts(i), ci_rmat)
174  else
175  do j = 1, scatci_input % opts(i) % num_target_sym
176  allocate(target_ci_hamiltonian::tgt_ci_hamiltonian)
177  call dispatchmatrixanddiagonalizer(scatci_input % opts(i) % diagonalizer_choice, &
178  scatci_input % opts(i) % force_serial, &
179  scatci_input % opts(i) % num_ci_target_sym(j), &
180  ci_rmat(j) % nstat, &
181  matrix_elements, &
182  diagonalizer, &
183  integrals, &
184  scatci_input % opts(i) % hamiltonian_unit)
185 
186  call matrix_elements % construct(target_hamiltonian + j)
187  call matrix_elements % set_options(scatci_input % opts(i))
188 
189  call tgt_ci_hamiltonian % construct(scatci_input % opts(i), csfs, orbitals, integrals)
190  call tgt_ci_hamiltonian % initialize(j)
191 
192  call master_timer % start_timer("Build CI Hamiltonian")
193  call tgt_ci_hamiltonian % build_hamiltonian(matrix_elements)
194  call master_timer % stop_timer("Build CI Hamiltonian")
195  call master_timer % report_timers
196 
197  call master_timer % start_timer("Diagonalization")
198  call diagonalizer % diagonalize(matrix_elements, ci_rmat(j) % nstat, ci_rmat(j), .true., &
199  scatci_input % opts(i), integrals)
200  call master_timer % stop_timer("Diagonalization")
201 
202  call ci_rmat(j) % print
203  call master_timer % report_timers
204 
205  call matrix_elements % destroy
206  if (process_grid % grank == master) close (scatci_input % opts(i) % hamiltonian_unit)
207  deallocate(matrix_elements, diagonalizer)
208  deallocate(tgt_ci_hamiltonian)
209  end do
210  end if
211  end if
212 
213  if (scatci_input % opts(i) % diagonalization_flag /= no_diagonalization) then
214  call dispatchmatrixanddiagonalizer(scatci_input % opts(i) % diagonalizer_choice, &
215  scatci_input % opts(i) % force_serial, &
216  scatci_input % opts(i) % contracted_mat_size, &
217  scatci_input % opts(i) % num_eigenpairs, &
218  matrix_elements, &
219  diagonalizer, &
220  integrals, &
221  scatci_input % opts(i) % hamiltonian_unit)
222  call matrix_elements % construct(main_hamiltonian)
223  call matrix_elements % set_options(scatci_input % opts(i))
224  else
225  call integrals % write_matrix_header(scatci_input % opts(i) % hamiltonian_unit, &
226  scatci_input % opts(i) % contracted_mat_size)
227  allocate(writermatrix::matrix_elements)
228  call matrix_elements % construct(main_hamiltonian)
229  call matrix_elements % set_options(scatci_input % opts(i))
230  end if
231 
232  call matrix_elements % exclude_row_column(scatci_input % opts(i) % exclude_rowcolumn)
233 
234  ! We are doing a scattering calculation
235  if (scatci_input % opts(i) % do_ci_contraction()) then
236  call enrgmx % construct(scatci_input % opts(i), csfs, orbitals, integrals)
237  call enrgmx % initialize(ci_rmat)
238  call master_timer % start_timer("C-Hamiltonian Build")
239  call enrgmx % build_hamiltonian(matrix_elements)
240  call master_timer % stop_timer("C-Hamiltonian Build")
241  deallocate (ci_rmat)
242  else
243  call enrgms % construct(scatci_input % opts(i), csfs, orbitals, integrals)
244  call master_timer % start_timer("Target-Hamiltonian Build")
245  call enrgms % build_hamiltonian(matrix_elements)
246  call master_timer % stop_timer("Target-Hamiltonian Build")
247  end if
248 
249  ! Lets deallocate all the csfs as well
250  deallocate(csfs)
251 
252  ! Clear the integrals, we dont need them anymore
253  call configuration_manager % finalize
254 
255  call master_memory % print_memory_report
256  write (stdout, "('Matrix N is ',i8)") matrix_elements % get_matrix_size()
257  write (stdout, "('Num elements is ',i14)") matrix_elements % get_size()
258  call master_timer % report_timers
259 
260  if (scatci_input % opts(i) % diagonalization_flag /= no_diagonalization) then
261  scatci_input % opts(i) % num_eigenpairs = min(scatci_input % opts(i) % num_eigenpairs, &
262  matrix_elements % get_matrix_size())
263  call solutions(i) % construct(scatci_input % opts(i))
264  call master_timer % start_timer("Diagonalization")
265  call diagonalizer % diagonalize(matrix_elements, &
266  scatci_input % opts(i) % num_eigenpairs, &
267  solutions(i), &
268  .false., &
269  scatci_input % opts(i), &
270  integrals)
271  call master_timer % stop_timer("Diagonalization")
272  else
273  write (stdout, '(/,"Diagonalization not selected")')
274  num_mat_elms = matrix_elements % get_size()
275 
276  call master_timer % start_timer("Matrix Save")
277  !call matrix_elements % save(SCATCI_input % hamiltonian_unit, SCATCI_input % num_matrix_elements_per_rec)
278  call master_timer % stop_timer("Matrix Save")
279 
280  if (process_grid % grank == master) call ham_pars_io_wrapper(scatci_input % opts(i), .true., num_mat_elms)
281 
282  write (stdout, '(/,"Parameters written to: ham_data")')
283  end if
284 
285  if (process_grid % grank == master) close (scatci_input % opts(i) % hamiltonian_unit)
286 
287  call matrix_elements % destroy
288  call master_timer % stop_timer("Wall Time")
289  call master_timer % report_timers
290 
291  deallocate (diagonalizer)
292  deallocate (matrix_elements)
293 
294  end do symmetry_loop
295 
296  ! use solutions in outer interface etc.
297  call postprocess(scatci_input, solutions)
298 
299  ! release (possibly shared) memory occupied by the integral storage
300  if (associated(integrals)) then
301  call integrals % finalize
302  nullify (integrals)
303  end if
304 
305  ! cleanly exit MPI
306  call mpi_mod_finalize
307 
308 end program mpi_scatci
contracted_hamiltonian_module::contracted_hamiltonian
Computation of Hamiltonian.
Definition: Contracted_Hamiltonian_module.f90:65
csf_module
CSF module.
Definition: CSF_module.F90:34
writermatrix_module
Write matrix module.
Definition: WriterMatrix_module.f90:31
memorymanager_module
Memory manager module.
Definition: MemoryManager_module.f90:31
postprocessing_module
Further processing of the diagonalization results.
Definition: Postprocessing_module.F90:32
options_module::ham_pars_io_wrapper
subroutine, public ham_pars_io_wrapper(option, write_ham_pars, nelms)
Read/write information about Hamiltonian.
Definition: Options_module.f90:1063
solutionhandler_module::solutionhandler
Solution writer.
Definition: SolutionHandler_module.f90:56
timing_module
Timing module.
Definition: Timing_Module.f90:28
ci_hamiltonian_module::target_ci_hamiltonian
Hamiltonian type.
Definition: CI_Hamiltonian_module.f90:56
target_rmat_ci_module::target_rmat_ci
This class handles the storage of target CI coefficients.
Definition: Target_RMat_CI_module.f90:52
consts_mpi_ci::no_ci_target
integer, parameter no_ci_target
No Ci target contraction (target only run)
Definition: consts_mpi_ci.f90:61
uncontracted_hamiltonian_module::uncontracted_hamiltonian
Definition: Uncontracted_Hamiltonian_module.f90:45
csf_module::csfobject
Single configuration state function.
Definition: CSF_module.F90:135
writermatrix_module::writermatrix
Matrix associated with a disk drive.
Definition: WriterMatrix_module.f90:56
ci_hamiltonian_module
Target CI Hamiltonian module.
Definition: CI_Hamiltonian_module.f90:31
parallelization_module::process_grid
type(processgrid) process_grid
Definition: Parallelization_module.F90:80
dispatcher_module::dispatchintegral
subroutine, public dispatchintegral(sym_group_flag, ukrmolplus_integrals, integral)
This subroutine dispatches the correct integral class based on simple parameters.
Definition: Dispatcher_module.F90:243
basematrix_module
Base matrix module.
Definition: BaseMatrix_module.f90:28
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
mpi_scatci
program mpi_scatci
MPI-SCATCI main program.
Definition: mpi_scatci.f90:41
contracted_hamiltonian_module
Contracted Hamiltonian module.
Definition: Contracted_Hamiltonian_module.f90:33
consts_mpi_ci
MPI SCATCI Constants module.
Definition: consts_mpi_ci.f90:31
options_module::optionsset
MPI-SCATCI input.
Definition: Options_module.f90:202
orbital_module
Target CI Hamiltonian module.
Definition: Orbital_module.f90:31
dispatcher_module
Dispatcher module.
Definition: Dispatcher_module.F90:32
basematrix_module::basematrix
Base matrix type.
Definition: BaseMatrix_module.f90:47
baseintegral_module
Base integral module.
Definition: BaseIntegral_module.f90:28
parallelization_module
Distribution of processes into a grid.
Definition: Parallelization_module.F90:29
dispatcher_module::dispatchmatrixanddiagonalizer
subroutine, public dispatchmatrixanddiagonalizer(diag_choice, force_serial, matrix_size, number_of_eigenpairs, matrix, diagonalizer, integral, matrix_io)
This subroutine determines which matrix format/diagonalizer pair to be used by SCATCI in build diagon...
Definition: Dispatcher_module.F90:113
csf_module::csforbital
Single determinant.
Definition: CSF_module.F90:118
diagonalizer_module
Base type for all diagonalizers.
Definition: Diagonalizer_module.f90:30
diagonalizer_module::basediagonalizer
Definition: Diagonalizer_module.f90:34
baseintegral_module::baseintegral
Definition: BaseIntegral_module.f90:41
uncontracted_hamiltonian_module
Uncontracted Hamiltonian module.
Definition: Uncontracted_Hamiltonian_module.f90:28
consts_mpi_ci::main_hamiltonian
integer, parameter main_hamiltonian
Hamiltonain is a scattering one.
Definition: consts_mpi_ci.f90:137
solutionhandler_module
Solution handler module.
Definition: SolutionHandler_module.f90:31
options_module
Options module.
Definition: Options_module.f90:41
timing_module::master_timer
type(timer), public master_timer
Definition: Timing_Module.f90:74
orbital_module::orbitaltable
This class generates the molecular and spin orbitals, stores them and generates symblic elements from...
Definition: Orbital_module.f90:69
csf_module::csfmanager
Configuration state function factory.
Definition: CSF_module.F90:63
consts_mpi_ci::no_diagonalization
integer, parameter no_diagonalization
No diagonalization.
Definition: consts_mpi_ci.f90:85
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
dispatcher_module::initialize_libraries
subroutine, public initialize_libraries
Initialize libraries used by MPI-SCATCI.
Definition: Dispatcher_module.F90:86
consts_mpi_ci::target_hamiltonian
integer, parameter target_hamiltonian
Hamiltonain is target one.
Definition: consts_mpi_ci.f90:139
postprocessing_module::postprocess
subroutine, public postprocess(SCATCI_input, solutions)
Full post-processing pass.
Definition: Postprocessing_module.F90:104