generate_fci_core Subroutine

public subroutine generate_fci_core(ilut_list, space_size)

Arguments

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

Contents

Source Code


Source Code

    subroutine generate_fci_core(ilut_list, space_size)

        ! In/Out: ilut_list - List of determinants generated.
        ! In/Out: space_size - Number of determinants in the generated space.
        !             If ilut_list is not empty on input and you want to keep
        !             the states already in it, then on input space_size should
        !             be equal to the number of states to be kept in ilut_list,
        !             and new states will be added in from space_size+1.
        !             Otherwise, space_size must equal 0 on input.
        !             On output space_size will equal the total number of
        !             generated plus what space_size was on input.


        integer(n_int), intent(inout) :: ilut_list(0:, :)
        integer, intent(inout) :: space_size

        integer, allocatable :: nI_list(:, :)
        integer(n_int) :: ilut(0:NIfTot)
        integer :: temp(1, 1), hf_ind, ndets, i
        character(*), parameter :: t_r = "generate_fci_core"

        if (.not. tSymSet) call stop_all(t_r, "To use the 'FCI-CORE' option you must also choose the symmetry sector of &
                                              &space to be generated by using the 'SYM' option.")

        if (.not. (tSpn .or. tGUGACore)) then
            call stop_all(t_r, "To use the FCI-CORE option you must fix the total spin on &
                               &z-axis, e.g., by using 'SPIN-RESTRICT' option for a SD &
                               &basis or 'GUGA' option for a GT basis.")
        end if

        ! Count the total number of determinants.
        call gndts(nel, nbasis, BRR, nBasisMax, temp, .true., G1, tSpn, lms, tParity, SymRestrict, ndets, hf_ind)
        allocate(nI_list(nel, ndets))
        if (size(ilut_list, 2) < ndets) then
            call stop_all(t_r, 'The ilut_list argument is too small. Probably SpawnedParts was allocated too small &
                    & Please increase memoryfacspawn and memoryfacpart')
        end if
        ! Generate and store all the determinants in nI_list.
        call gndts(nel, nbasis, BRR, nBasisMax, nI_list, .false., G1, tSpn, lms, tParity, SymRestrict, ndets, hf_ind)

        if (tGUGACore) then
            do i = 1, ndets
                call EncodeBitDet(nI_list(:, i), ilut)
                if (isProperCSF_flexible(ilut, STOT, nel)) call add_state_to_space(ilut, ilut_list, space_size, nI_list(:, i))
            end do
        else
            do i = 1, ndets
                call EncodeBitDet(nI_list(:, i), ilut)
                call add_state_to_space(ilut, ilut_list, space_size, nI_list(:, i))
            end do
        end if

    end subroutine generate_fci_core