subroutine InitHamiltonianCalc(this, print_info, hamil_type, max_subspace_size, t_store_subspace_basis, t_orthogonalise)
use direct_ci, only: create_ham_diag_direct_ci
use FciMCData, only: davidson_ras, davidson_classes, davidson_strings
use ras, only: find_ras_size
use util_mod, only: int_fmt
type(HamiltonianCalcType), intent(inout) :: this
logical, intent(in) :: print_info
integer :: i, mem_reqd, residual_mem_reqd, ierr, space_size
integer, intent(in) :: hamil_type, max_subspace_size
logical, intent(in) :: t_store_subspace_basis, t_orthogonalise
character(*), parameter :: t_r = "InitHamiltonianCalc"
integer(MPIArg) :: mpi_temp
real(dp), allocatable :: hamil_diag_temp(:)
this%hamil_type = hamil_type
this%t_store_subspace_basis = t_store_subspace_basis
this%t_orthogonalise = t_orthogonalise
associate ( &
space_size => this%space_size, &
hfindex => this%hfindex &
)
! allocate and define the hamiltonian diagonal, if not done so already.
if (.not. allocated(hamil_diag)) then
if (hamil_type == direct_ci_type) then
allocate(hamil_diag(space_size), stat=ierr)
call logmemalloc("hamil_diag", space_size, 8, t_r, hdiagtag, ierr)
call create_ham_diag_direct_ci(davidson_ras, davidson_classes, davidson_strings, hamil_diag)
else if (hamil_type == full_hamil_type) then
space_size = size(hamiltonian, 1)
allocate(hamil_diag(space_size), stat=ierr)
call logmemalloc("hamil_diag", space_size, 8, t_r, hdiagtag, ierr)
do i = 1, space_size
hamil_diag(i) = hamiltonian(i, i)
end do
else if (hamil_type == parallel_sparse_hamil_type .or. hamil_type == sparse_hamil_type) then
! in the case of sparse implementations, the diagonal should be
! created when the hamiltonian itself is.
call stop_all("t_r", "the diagonal of the hamiltonian has not been allocated. &
&cannot perform Davidson or Lanczos calculation.")
end if
end if
space_size = size(hamil_diag)
! may not exceed size of space)
this%max_subspace_size = max_subspace_size!min(max_subspace_size, space_size)
if (hamil_type == parallel_sparse_hamil_type) then
allocate(this%space_sizes(0:nProcessors - 1))
allocate(this%partial_H_ket_disps(0:nProcessors - 1))
allocate(this%partial_H_ket(space_size))
mpi_temp = int(space_size, MPIArg)
call MPIAllGather(mpi_temp, this%space_sizes, ierr)
! The total space size across all processors.
space_size = int(sum(this%space_sizes))
allocate(hamil_diag_temp(space_size))
this%partial_H_ket_disps(0) = 0
do i = 1, nProcessors - 1
this%partial_H_ket_disps(i) = this%partial_H_ket_disps(i - 1) &
+ this%space_sizes(i - 1)
end do
call MPIGatherV(hamil_diag, hamil_diag_temp, this%space_sizes, this%partial_H_ket_disps, ierr)
if (iProcIndex == root) then
safe_realloc(hamil_diag, (space_size))
hamil_diag = hamil_diag_temp
end if
deallocate(hamil_diag_temp)
end if
if (print_info) then
write(stdout, '(1x,"number of determinants in total:",'//int_fmt(space_size, 1)//')') space_size; call neci_flush(stdout)
end if
if (tCalcHFIndex) then
hfindex = maxloc((-hamil_diag), 1)
end if
! the memory required to allocate each of basis_vectors and
! multipied_basis_vectors, in mb.
mem_reqd = max_subspace_size * this%space_size * 8 / 1000000
! allocate the necessary arrays:
if (t_store_subspace_basis) then
safe_calloc(this%basis_vectors, (this%space_size, max_subspace_size), 0.0_dp)
if (print_info) then
write(stdout, '(1x,"allocating array to hold subspace vectors (",'//int_fmt(mem_reqd, 0)//',1x,"mb).")') mem_reqd
call neci_flush(stdout)
end if
end if
safe_calloc(this%projected_hamil, (max_subspace_size, max_subspace_size), 0.0_dp)
safe_calloc(this%projected_hamil_work, (max_subspace_size, max_subspace_size), 0.0_dp)
if (iprocindex /= root) then
safe_malloc(this%temp_in, (this%space_size))
safe_malloc(this%temp_out, (this%space_size))
end if
this%skip_calc = .false.
if (print_info) then
write(stdout, '(1x,"Hamiltonian calculation setup complete.",/)')
end if
call neci_flush(stdout)
end associate
end subroutine InitHamiltonianCalc