An API interfacing function for generate_excitation to the rest of NECI:
Requires guga_bitRepOps::current_csf_i to be set according to the ilutI.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | nI(nel) | |||
integer(kind=n_int), | intent(in) | :: | ilutI(0:NIfTot) | |||
integer, | intent(out) | :: | nJ(nel) | |||
integer(kind=n_int), | intent(out) | :: | ilutJ(0:NifTot) | |||
integer, | intent(in) | :: | exFlag | |||
integer, | intent(out) | :: | ic | |||
integer, | intent(out) | :: | ex(2,maxExcit) | |||
logical, | intent(out) | :: | tParity | |||
real(kind=dp), | intent(out) | :: | pGen | |||
real(kind=dp), | intent(out) | :: | hel | |||
type(excit_gen_store_type), | intent(inout), | target | :: | store | ||
integer, | intent(in), | optional | :: | run |
subroutine gen_excit_k_space_hub(nI, ilutI, nJ, ilutJ, exFlag, ic, &
ex, tParity, pGen, hel, store, run)
!! An API interfacing function for generate_excitation to the rest of NECI:
!!
!! Requires guga_bitRepOps::current_csf_i to be set according to the ilutI.
integer, intent(in) :: nI(nel), exFlag
integer(n_int), intent(in) :: ilutI(0:NIfTot)
integer, intent(out) :: nJ(nel), ic, ex(2, maxExcit)
integer(n_int), intent(out) :: ilutJ(0:NifTot)
real(dp), intent(out) :: pGen
logical, intent(out) :: tParity
HElement_t(dp), intent(out) :: hel
type(excit_gen_store_type), intent(inout), target :: store
integer, intent(in), optional :: run
#ifdef DEBUG_
character(*), parameter :: this_routine = "gen_excit_k_space_hub"
#endif
real(dp) :: p_elec, p_orb
integer :: elecs(2), orbs(2), src(2)
type(ExcitationInformation_t) :: excitInfo
integer(n_int) :: ilutGj(0:nifguga)
unused_var(exFlag); unused_var(store); unused_var(run)
hel = h_cast(0.0_dp)
ic = 0
pgen = 0.0_dp
! i first have to choose an electron pair (ij) at random
! but with the condition that they have to have opposite spin!
call pick_spin_opp_elecs(nI, elecs, p_elec)
src = nI(elecs)
call pick_ab_orbitals_hubbard(nI, ilutI, src, orbs, p_orb)
if (orbs(1) == ABORT_EXCITATION) then
nJ(1) = ABORT_EXCITATION
pgen = 0.0_dp
return
end if
ic = 2
! and make the excitation
call make_double(nI, nJ, elecs(1), elecs(2), orbs(1), orbs(2), ex, tParity)
ilutJ = make_ilutJ(ilutI, ex, 2)
! i think in both the electrons and the orbitals i have twice the
! probability to pick them
! already modified in the orbital picker..
pgen = p_elec * p_orb
! try implementing the crude guga excitation approximation via the
! determinant excitation generator
if (tGen_guga_crude) then
call convert_ilut_toGUGA(ilutJ, ilutGj)
if (.not. isProperCSF_ilut(ilutGJ, .true.)) then
nJ(1) = 0
pgen = 0.0_dp
return
end if
ASSERT(is_compatible(ilutI, current_csf_i))
call calc_guga_matrix_element(ilutI, current_csf_i, ilutJ, CSF_Info_t(ilutJ), excitInfo, hel, .true.)
if (abs(hel) < EPS) then
nJ(1) = 0
pgen = 0.0_dp
end if
global_excitinfo = excitInfo
return
end if
#ifdef DEBUG_
if (.not. isvaliddet(nI, nel)) then
if (nJ(1) /= 0) then
print *, "nI: ", nI
print *, "nJ: ", nJ
print *, "src: ", ex(1, :)
print *, "tgt: ", ex(2, :)
end if
end if
#endif
end subroutine gen_excit_k_space_hub