generate_trip_determinants Subroutine

private subroutine generate_trip_determinants(ilut_list, space_size, only_keep_conn)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(inout) :: ilut_list(0:,:)
integer, intent(inout) :: space_size
logical, intent(in) :: only_keep_conn

Contents


Source Code

    subroutine generate_trip_determinants(ilut_list, space_size, only_keep_conn)
        ! Generate a list of all singles, doubles and triples
        integer(n_int), intent(inout) :: ilut_list(0:, :)
        integer, intent(inout) :: space_size
        logical, intent(in) :: only_keep_conn

        integer :: i, j, k
        integer :: a, b, c
        integer :: ex(2, 3), nI(nel)
        integer(n_int) :: ilut_ex(0:NifTot)
        HElement_t(dp) :: HEl
        ! add all triples of the HF
        ! enumerate them as follows: take three electrons i,j,k
        ! with i>j>k
        ! then, pick three (unocc) orbitals a,b,c with a>b>c
        ! => triple excitation
        do i = 1, nel
            ex(1, 1) = HFDet(i)
            do j = 1, i - 1
                ex(1, 2) = HFDet(j)
                do k = 1, j - 1
                    ex(1, 3) = HFDet(k)
                    do a = 1, nBasis
                        ! for the target orbs, only take unoccupied
                        if (IsNotOcc(ilutHF, a)) then
                            ex(2, 3) = a
                            do b = 1, a - 1
                                if (IsNotOcc(ilutHF, b)) then
                                    ex(2, 2) = b
                                    do c = 1, b - 1
                                        if (IsNotOcc(ilutHF, c)) then
                                            ex(2, 1) = c
                                            ! create the triple excitation with these elecs/orbs
                                            ilut_ex = make_ilutJ(ilutHF, ex, 3)
                                            ! if enabled, only keep connected determinants
                                            if (only_keep_conn) then
                                                HEl = get_helement(HFDet, nI, ilutHF, ilut_ex)
                                                if (abs(HEl) < eps) cycle
                                            end if

                                            ! definitely keep only determinants in the same symmetry-sector
                                            if (.not. IsSymAllowedExcitMat(ex, 3)) cycle
                                            call decode_bit_det(nI, ilut_ex)
                                            call add_state_to_space(ilut_ex, ilut_list, space_size, nI)
                                        end if
                                    end do
                                end if
                            end do
                        end if
                    end do
                end do
            end do
        end do

    end subroutine generate_trip_determinants