fill_csf_i Subroutine

public pure subroutine fill_csf_i(ilut, csf_i)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut(0:GugaBits%len_tot)
type(CSF_Info_t), intent(inout) :: csf_i

Contents

Source Code


Source Code

    pure subroutine fill_csf_i(ilut, csf_i)
        ! routine which sets up all the additional csf information, like
        ! stepvector, b vector, occupation etc. in various formats in one
        ! place
        ! and combine all the necessary calcs. into one loop instead of
        ! the seperate ones..
        integer(n_int), intent(in) :: ilut(0:GugaBits%len_tot)
        type(CSF_Info_t), intent(inout) :: csf_i
        debug_function_name("fill_csf_i")

        integer :: i, step, b_int
        real(dp) :: b_real, cum_sum

        ASSERT(isProperCSF_ilut(ilut))
        ASSERT(allocated(csf_i%stepvector))
        ASSERT(allocated(csf_i%B_real))
        ASSERT(allocated(csf_i%Occ_real))
        ASSERT(allocated(csf_i%B_int))
        ASSERT(allocated(csf_i%Occ_int))

        csf_i%stepvector = 0
        csf_i%B_real = 0.0_dp
        csf_i%Occ_real = 0.0_dp
        csf_i%B_int = 0
        csf_i%Occ_int = 0

        b_real = 0.0_dp
        b_int = 0


        ! TODO(@Oskar): Use these functions instead
        ! currentB_ilut = calcB_vector_ilut(ilut)
        ! currentOcc_ilut = calcOcc_vector_ilut(ilut)
        ! currentOcc_int = calcOcc_vector_int(ilut)
        ! current_stepvector = calcStepVector(ilut)
        ! currentB_int = calcB_vector_int(ilut)

        ! also create a fake cum-list of the non-doubly occupied orbitals
        csf_i%cum_list = 0.0_dp
        cum_sum = 0.0_dp

        do i = 1, nSpatOrbs

            step = getStepvalue(ilut, i)

            csf_i%stepvector(i) = step

            select case (step)

            case (0)

                csf_i%Occ_real(i) = 0.0_dp
                csf_i%Occ_int(i) = 0

                cum_sum = cum_sum + 1.0_dp

            case (1)

                csf_i%Occ_real(i) = 1.0_dp
                csf_i%Occ_int(i) = 1

                b_real = b_real + 1.0_dp
                b_int = b_int + 1

                cum_sum = cum_sum + 1.0_dp

            case (2)

                csf_i%Occ_real(i) = 1.0_dp
                csf_i%Occ_int(i) = 1

                b_real = b_real - 1.0_dp
                b_int = b_int - 1

                cum_sum = cum_sum + 1.0_dp

            case (3)

                csf_i%Occ_real(i) = 2.0_dp
                csf_i%Occ_int(i) = 2

            end select

            csf_i%B_real(i) = b_real
            csf_i%B_int(i) = b_int

            csf_i%cum_list(i) = cum_sum
        end do

    end subroutine fill_csf_i