frozen_single_entry Function

private function frozen_single_entry(indices, prefactor) result(orbs)

For a given set of 6 orbital indices with two pairs of frozen orbitals, returns the orbital indices of the corresponding single excitation and the prefactor for the matrix element due to spin @param[in] indices array of size 6 with orbital indices, four of which have to be repeated frozen indices @param[out] prefactor on return, the prefactor of the matrix element when added to the one-body terms @return orbs the two non-frozen orbitals

Arguments

Type IntentOptional Attributes Name
integer(kind=int64), intent(in) :: indices(num_inds)
real(kind=dp), intent(out) :: prefactor

Return Value integer, (2)


Contents

Source Code


Source Code

    function frozen_single_entry(indices, prefactor) result(orbs)
        integer(int64), intent(in) :: indices(num_inds)
        real(dp), intent(out) :: prefactor

        integer :: orbs(2)
        integer :: f_orbs(num_inds - 2), ct
        integer :: uf_one, uf_two, directs
        logical :: unfrozen(num_inds)

        orbs = 0
        unfrozen = indices > 0
        f_orbs = int(pack(indices,.not. unfrozen))
        ! For the single entry, there are only two unfrozen orbs, find these and return them
        orbs = int(pack(indices, unfrozen))

        ! Again, two things to check: The number of direct orbitals and the permutation

        ! At this point, there is no unpaired frozen orbital, so there are two options:
        ! a) there are two different orbs
        prefactor = 1.0_dp
        if(count(f_orbs(1) == f_orbs) < size(f_orbs)) then
            ! We have to count how many potential spin configurations contribute

            ! The parity is counting the number of permutations required to end up
            ! with a completely direct expression
            ! -> The parity is odd if and only if the number of direct excitations is exactly one
            ! (then, one permutation is required to fix the two non-direct excits, else,
            ! (either zero or two are needed)
            ! -> count every direct excitation with a factor -1
            directs = 0
            ! Count the number of direct excitations
            do ct = 1, step
                ! Each direct repeated orbital doubles the number of spin configs
                ! (the case of four identical orbs is excluded)
                if(is_repeated_pair(indices, ct) .or. (unfrozen(ct) .and. unfrozen(ct + step))) then
                    directs = directs + 1
                endif
            end do
            select case(directs)
                ! If there are none, there is no freedom for spin choice, but the LMat entry
                ! appears twice in the matrix element evaluation due to symmetry if the
                ! unfrozen entry is repeated
            case(0)
                if(orbs(1) == orbs(2)) prefactor = 2.0_dp
                ! If there is one, two possible spin configs contribute, and we have odd parity
            case(1)
                prefactor = -2.0_dp
            case(3)
                ! If there are three (just two is not possible), there are two spin degrees of
                ! freedom -> factor of 4
                prefactor = 4.0_dp
            end select
        else
            ! b) All frozen orbs are the same -> only one spin config is allowed, only parity of
            ! the permutation counts
            ! Here, only a direct exctiation has even parity (same rule as above)

            ! Get the posisitons of the unfrozen orbs
            uf_one = custom_findloc(unfrozen, .true., back=.false.)
            uf_two = custom_findloc(unfrozen, .true., back=.true.)

            if(.not. is_direct(uf_one, uf_two)) prefactor = -1.0_dp * prefactor
        end if

    end function frozen_single_entry