subroutine guess_huckel(infos)
use precision, only: dp
use types, only: information
use io_constants, only: IW
use oqp_tagarray_driver
use basis_tools, only: basis_set
use guess, only: get_ab_initio_density
use mathlib, only: matrix_invsqrt
use huckel, only: huckel_guess
use util, only: measure_time
use messages, only: show_message, WITH_ABORT
use strings, only: Cstring, fstring
use printing, only: print_module_info
use oqp_tagarray_driver
use iso_c_binding, only: c_char
use parallel, only: par_env_t
implicit none
character(len=*), parameter :: subroutine_name = "guess_huckel"
type(information), target, intent(inout) :: infos
integer :: i, nbf, nbf2
type(basis_set), pointer :: basis
type(basis_set) :: huckel_basis
character(len=:), allocatable :: basis_file
logical :: err
integer , parameter :: root = 0
type(par_env_t) :: pe
! tagarray
real(kind=dp), contiguous, pointer :: &
Smat(:), &
dmat_a(:), mo_a(:,:), mo_energy_a(:), &
dmat_b(:), mo_b(:,:), mo_energy_b(:)
character(len=*), parameter :: tags_alpha(3) = (/ character(len=80) :: &
OQP_DM_A, OQP_E_MO_A, OQP_VEC_MO_A /)
character(len=*), parameter :: tags_beta(3) = (/ character(len=80) :: &
OQP_DM_B, OQP_E_MO_B, OQP_VEC_MO_B /)
character(len=*), parameter :: tags_general(2) = (/ character(len=80) :: &
OQP_SM, OQP_hbasis_filename /)
character(len=1,kind=c_char), contiguous, pointer :: basis_filename(:)
!
! Section of Tagarray for the log filename
! We are getting lot file name from Python via tagarray
!
call data_has_tags(infos%dat, tags_general, module_name, subroutine_name, with_abort)
call tagarray_get_data(infos%dat, OQP_hbasis_filename, basis_filename)
allocate(character(ubound(basis_filename,1)) :: basis_file)
do i = 1, ubound(basis_filename,1)
basis_file(i:i) = basis_filename(i)
end do
!
! Files open
! 1. XYZ: Read : Geometric data, ATOMS
! 3. LOG: Read Write: Main output file
!
open (unit=IW, file=infos%log_filename, position="append")
!
call print_module_info('Guess_Huckel','Initial guess using Huckel theory')
!
! Readings
! load basis set
basis => infos%basis
call pe%init(infos%mpiinfo%comm, infos%mpiinfo%usempi)
if (pe%rank == root) then
call huckel_basis%from_file(basis_file, infos%atoms, err)
end if
! Checking error of basis set reading..
infos%control%basis_set_issue = err
call pe%bcast(infos%control%basis_set_issue, 1)
basis%atoms => infos%atoms
! Allocate H, S ,T and D matrices
nbf = basis%nbf
nbf2 =nbf*(nbf+1)/2
! clean data
call infos%dat%remove_records(tags_alpha)
call infos%dat%remove_records(tags_beta)
! load general data
call data_has_tags(infos%dat, tags_general, module_name, subroutine_name, WITH_ABORT)
call tagarray_get_data(infos%dat, OQP_SM, smat)
! allocate alpha
call infos%dat%reserve_data(OQP_DM_A, TA_TYPE_REAL64, nbf2, comment=OQP_DM_A_comment)
call infos%dat%reserve_data(OQP_E_MO_A, TA_TYPE_REAL64, nbf, comment=OQP_E_MO_A_comment)
call infos%dat%reserve_data(OQP_VEC_MO_A, TA_TYPE_REAL64, nbf*nbf, (/ nbf, nbf /), comment=OQP_VEC_MO_A_comment)
! load alpha data
call data_has_tags(infos%dat, tags_alpha, module_name, subroutine_name, WITH_ABORT)
call tagarray_get_data(infos%dat, OQP_DM_A, dmat_a)
call tagarray_get_data(infos%dat, OQP_E_MO_A, mo_energy_a)
call tagarray_get_data(infos%dat, OQP_VEC_MO_A, mo_a)
! UHF/ROHF
if (infos%control%scftype >= 2) then
! allocate beta
call infos%dat%reserve_data(OQP_DM_B, TA_TYPE_REAL64, nbf2, comment=OQP_DM_B_comment)
call infos%dat%reserve_data(OQP_E_MO_B, TA_TYPE_REAL64, nbf, comment=OQP_E_MO_B_comment)
call infos%dat%reserve_data(OQP_VEC_MO_B, TA_TYPE_REAL64, nbf*nbf, (/ nbf, nbf /), comment=OQP_VEC_MO_B_comment)
! load beta
call data_has_tags(infos%dat, tags_beta, module_name, subroutine_name, WITH_ABORT)
call tagarray_get_data(infos%dat, OQP_DM_B, dmat_b)
call tagarray_get_data(infos%dat, OQP_E_MO_B, mo_energy_b)
call tagarray_get_data(infos%dat, OQP_VEC_MO_B, mo_b)
end if
! Calculate Huckel MO
if (pe%rank == root) then
call huckel_guess(Smat, MO_A, infos, basis, huckel_basis)
endif
! For ROHF/UHF
if (INFOS%control%scftype >= 2) MO_B = MO_A
!
! Calculate Density Matrix
if (pe%rank == root) then
! RHF
if (infos%control%scftype == 1) then
call get_ab_initio_density(Dmat_A, MO_A, infos=infos, basis=basis)
! ROHF/UHF
else
call get_ab_initio_density(Dmat_A, MO_A, Dmat_B, MO_B, infos, basis)
endif
endif
! Broadcast MO and density matrices to all processes
call pe%bcast(MO_A, nbf*nbf)
if (infos%control%scftype >= 2) then
call pe%bcast(MO_B, nbf*nbf)
endif
! Broadcast the density matrices to all processes
if (infos%control%scftype == 1) then
call pe%bcast(Dmat_A, nbf2)
else
call pe%bcast(Dmat_A, nbf2)
call pe%bcast(Dmat_B, nbf2)
endif
call pe%barrier()
write (iw, '(/x,a,/)') '...... End of initial orbital guess ......'
call measure_time(print_total=1, log_unit=iw)
close(iw)
end subroutine guess_huckel