calc_guga_matrix_element Subroutine

public pure subroutine calc_guga_matrix_element(ilutI, csf_i, ilutJ, csf_j, excitInfo, mat_ele, t_hamil, rdm_ind, rdm_mat)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilutI(0:niftot)
type(CSF_Info_t), intent(in) :: csf_i
integer(kind=n_int), intent(in) :: ilutJ(0:niftot)
type(CSF_Info_t), intent(in) :: csf_j
type(ExcitationInformation_t), intent(out) :: excitInfo
real(kind=dp), intent(out) :: mat_ele
logical, intent(in) :: t_hamil
integer(kind=int_rdm), intent(out), optional, allocatable :: rdm_ind(:)
real(kind=dp), intent(out), optional, allocatable :: rdm_mat(:)

Contents


Source Code

    pure subroutine calc_guga_matrix_element(ilutI, csf_i, ilutJ, csf_j, excitInfo, mat_ele, t_hamil, &
                                        rdm_ind, rdm_mat)
        ! function which, given the 2 CSFs ilutI/J and the excitation
        ! information, connecting those 2, calculates the Hamiltionian
        ! matrix element between those 2
        ! use a flag to distinguish between only guga-mat_ele calculation
        ! and full hamiltonian matrix element calculation

        integer(n_int), intent(in) :: ilutI(0:niftot), ilutJ(0:niftot)
        type(CSF_Info_t), intent(in) :: csf_i, csf_j
        type(ExcitationInformation_t), intent(out) :: excitInfo
        HElement_t(dp), intent(out) :: mat_ele
        logical, intent(in) :: t_hamil
        integer(int_rdm), intent(out), allocatable, optional :: rdm_ind(:)
        real(dp), intent(out), allocatable, optional :: rdm_mat(:)
        character(*), parameter :: this_routine = "calc_guga_matrix_element"

        integer(n_int) :: tmp_i(0:nifguga)

        mat_ele = h_cast(0.0_dp)

        ASSERT(present(rdm_ind) .eqv. present(rdm_mat))

        ! check diagonal case first
        if (DetBitEQ(ilutI, ilutJ)) then

            call convert_ilut_toGUGA(ilutI, tmp_i)

            ! i think I should do a change here for the Heisenberg or
            ! tJ model..
            mat_ele = calcDiagMatEleGUGA_ilut(tmp_i)
            return
        end if

        excitInfo = identify_excitation(ilutI, ilutJ)

        ! more than a double excitation! leave with 0 matrix element
        if (.not. excitInfo%valid) return

        ! for the hubbard model implementation, depending if it is in the
        ! momentum- or real-space i can get out of here if we identify
        ! excitation types which are definetly 0
        ! i could do that more efficiently if we identify it already
        ! earlier but for now do it here!
        if (t_hamil) then
            ! make this check only if we want the hamiltonian matrix
            ! element. for general coupling coefficients (eg. for RDMs)
            ! i do need this contributions anyway
            if (t_new_hubbard) then
                if (treal .or. t_new_real_space_hubbard) then
                    ! only singles in the real-space hubbard!
                    if (excitInfo%typ /= excit_type%single) return
                else
                    ! only double excitations in the momentum-space hubbard!
                    if (excitInfo%typ == excit_type%single) return
                end if
            end if

            ! make the adjustment for the Heisenberg model
            if (t_heisenberg_model &
                .and. excitInfo%typ /= excit_type%fullstart_stop_mixed) then
                return
            end if

            if (t_tJ_model .and. &
                (.not. (excitInfo%typ == excit_type%single &
                 .or. excitInfo%typ == excit_type%fullstart_stop_mixed))) &
                return
        end if

        ! i think in the excitation identification i can not find out if the
        ! delta B value is abs>2 so i have to do that here.. or specific for
        ! the type of excitations below.. for singles its not allowed
        ! abs > 1 ..
        ! but i think for the double excitations i cannot do that generally
        ! since i have to check for the overlap and non-overlap regions
        ! specifically

        ! then i need a select case to specifically calculate all the
        ! different types of excitations
        ! essentially i can just mimick the stochastic excitation creation
        ! routines but with a fixed chosen excitation
        select case (excitInfo%typ)

        case (excit_type%single)
            ! pure single excitation

            ! but here i have to calculate all the double excitation
            ! influences which can lead to the same excitation(weights etc.)
            call calc_single_excitation_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                           t_hamil, rdm_ind, rdm_mat)

        case (excit_type%single_overlap_L_to_R)
            ! single overlap lowering into raising

            ! maybe i have to check special conditions on the overlap site.
            call calc_single_overlap_mixed_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                              t_hamil, rdm_ind, rdm_mat)

        case (excit_type%single_overlap_R_to_L)
            ! single overlap raising into lowering

            ! maybe i have to check special conditions on the overlap site.
            call calc_single_overlap_mixed_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                              t_hamil, rdm_ind, rdm_mat)

        case (excit_type%double_lowering)
            ! normal double lowering

            ! question is can i combine more functions here since i know
            ! both CSFs.. i think so!

            ! deal with order parameter for switched indices
            call calc_normal_double_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                       t_hamil, rdm_ind, rdm_mat)

        case (excit_type%double_raising)
            ! normal double raising

            ! here i have to deal with the order parameter for switched
            ! indices ..
            call calc_normal_double_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                       t_hamil, rdm_ind, rdm_mat)

        case (excit_type%double_L_to_R_to_L)
            ! lowering into raising into lowering

            ! can i combine these 4 similar excitations in one routine?
            ! deal with non-overlap if no spin-coupling changes!
            call calc_normal_double_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                       t_hamil, rdm_ind, rdm_mat)

        case (excit_type%double_R_to_L_to_R)
            ! raising into lowering into raising

            ! here i have to consider the non-overlap contribution if no
            ! spin-coupling changes in the overlap range
            call calc_normal_double_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                       t_hamil, rdm_ind, rdm_mat)

        case (excit_type%double_L_to_R)
            ! lowering into raising double

            ! consider non-overlap if no spin-coupling changes!
            call calc_normal_double_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                       t_hamil, rdm_ind, rdm_mat)

        case (excit_type%double_R_to_L)
            ! raising into lowering double

            ! here i also have to consider the non-overlap contribution if no
            ! spin-coupling changes in the overlap range
            call calc_normal_double_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                       t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstop_lowering)
            ! full-stop 2 lowering

            ! here only x0 matrix element in overlap range!
            ! also combine fullstop-alike
            call calc_fullstop_alike_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                        t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstop_raising)
            ! full-stop 2 raising

            ! here only x0 matrix elment in overlap range!
            call calc_fullstop_alike_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                        t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstop_L_to_R)
            ! full-stop lowering into raising

            ! here i have to consider all the singly occupied orbital
            ! influences ABOVE the last spin-coupling change
            ! this is more of a pain.. do later!
            ! finished the "easy" ones.. now to the annoying..
            call calc_fullstop_mixed_ex(ilutI, csf_i, ilutJ, csf_j, excitInfo, mat_ele, &
                                        t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstop_R_to_L)
            ! full-stop raising into lowering


            ! here i have to consider all the singly occupied orbital
            ! influences ABOVE the last spin-coupling change
            call calc_fullstop_mixed_ex(ilutI, csf_i, ilutJ, csf_j, excitInfo, mat_ele, &
                                        t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstart_lowering)
            ! full-start 2 lowering

            ! here only x0 matrix element in overlap range!
            call calc_fullstart_alike_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                         t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstart_raising)
            ! full-start 2 raising

            ! here only the x0-matrix in the overlap range (this implies no
            ! spin-coupling changes, but i already dealt with that! (hopefully!))
            call calc_fullstart_alike_ex(csf_i, csf_j, excitInfo, mat_ele, &
                                         t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullStart_L_to_R)
            ! full-start lowering into raising

            ! here i have to consider all the other singly occupied orbital
            ! influences BELOW the first spin-coupling change
            call calc_fullstart_mixed_ex(ilutI, csf_i, ilutJ, csf_j, excitInfo, mat_ele, &
                                         t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstart_R_to_L)
            ! full-start raising into lowering

            ! here i have to consider all the other singly occupied orbital
            ! influences BELOW the first spin-coupling change
            call calc_fullstart_mixed_ex(ilutI, csf_i, ilutJ, csf_j, excitInfo, mat_ele, &
                                         t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstart_stop_alike)
            ! full-start into full-stop alike

            ! here no spin-coupling changes are allowed!
            call calc_fullstart_fullstop_alike_ex(csf_i, excitInfo, &
                                                  mat_ele, t_hamil, rdm_ind, rdm_mat)

        case (excit_type%fullstart_stop_mixed)
            ! full-start into full-stop mixed

            ! here i have to consider all the singly occupied orbitals
            ! below the first spin-change and above the last spin change
            call calc_fullstart_fullstop_mixed_ex(ilutI, csf_i, ilutJ, csf_j, excitInfo, &
                                                  mat_ele, t_hamil, rdm_ind, rdm_mat)

        case default
            call stop_all(this_routine, "unexpected excitation type")

        end select

        if (tStoquastize) mat_ele = -abs(mat_ele)

        if (t_matele_cutoff .and. abs(mat_ele) < matele_cutoff) mat_ele = h_cast(0.0_dp)

    end subroutine calc_guga_matrix_element