singleStochasticEnd Subroutine

public subroutine singleStochasticEnd(csf_i, excitInfo, t)

Arguments

Type IntentOptional Attributes Name
type(CSF_Info_t), intent(in) :: csf_i
type(ExcitationInformation_t), intent(in) :: excitInfo
integer(kind=n_int), intent(inout) :: t(0:nifguga)

Contents

Source Code


Source Code

    subroutine singleStochasticEnd(csf_i, excitInfo, t)
        ! routine to end a stochastic excitation for a single generator
        type(CSF_Info_t), intent(in) :: csf_i
        type(ExcitationInformation_t), intent(in) :: excitInfo
        integer(n_int), intent(inout) :: t(0:nifguga)
        character(*), parameter :: this_routine = "singleStochasticEnd"

        integer :: ende, gen, deltaB
        real(dp) :: bVal, tempWeight

        ende = excitInfo%fullEnd
        deltaB = getDeltaB(t)
        tempWeight = extract_matrix_element(t, 1)
        bVal = csf_i%B_real(ende)
        gen = excitInfo%currentGen
        ! should be really similar to full excitation routine, except all the
        ! clean up stuff
        select case (csf_i%stepvector(ende))
        case (0)
            ! if it is zero it implies i have a lowering generator, otherwise
            ! i wouldnt even be here.
            if (deltaB == -1) then
                ! set alpha bit 0 -> 2
                set_orb(t, 2 * ende)

                tempWeight = getSingleMatrixElement(2, 0, deltaB, gen, bVal)

            else
                ! set beta bit 0 -> 1
                set_orb(t, 2 * ende - 1)

                tempWeight = getSingleMatrixElement(1, 0, deltaB, gen, bVal)

            end if

        case (3)
            ! if d = 3, it implies a raising generator or otherwise we
            ! would not be here ..
            if (deltaB == -1) then
                ! clear beta bit
                clr_orb(t, 2 * ende - 1)

                tempWeight = getSingleMatrixElement(2, 3, deltaB, gen, bVal)

            else
                ! clear alpha bit
                clr_orb(t, 2 * ende)

                tempWeight = getSingleMatrixElement(1, 3, deltaB, gen, bVal)

            end if

        case (1)

            if (t_approx_exchange .or. (t_approx_exchange_noninits .and. (.not. is_init_guga))) then
                if (deltaB /= -1) then
                    t = 0
                    return
                end if
            end if

            ASSERT(deltaB == -1)

            if (gen == gen_type%R) then
                clr_orb(t, 2 * ende - 1)

                tempWeight = getSingleMatrixElement(0, 1, deltaB, gen, bVal)

            else ! lowering generator

                ! for lowering gen: 1 -> 3: set alpha bit
                set_orb(t, 2 * ende)

                tempWeight = getSingleMatrixElement(3, 1, deltaB, gen, bVal)

            end if

        case (2)
            if (t_approx_exchange .or. (t_approx_exchange_noninits .and. (.not. is_init_guga))) then
                if (deltaB /= 1) then
                    t = 0
                    return
                end if
            end if
            ASSERT(deltaB == 1)

            if (gen == gen_type%R) then
                ! for raising: 2 -> 0: clear alpha bit
                clr_orb(t, 2 * ende)

                tempWeight = getSingleMatrixElement(0, 2, deltaB, gen, bVal)

            else ! lowering gen
                ! for lowering gen: 2 -> 3 : set beta bit
                set_orb(t, 2 * ende - 1)

                tempWeight = getSingleMatrixElement(3, 2, deltaB, gen, bVal)
            end if
        end select

        call update_matrix_element(t, tempWeight, 1)
        call update_matrix_element(t, tempWeight, 2)

        ! do excitaiton abortion
        if (near_zero(extract_matrix_element(t, 1)) .and. &
            near_zero(extract_matrix_element(t, 2))) then
            t = 0_n_int
        end if

    end subroutine singleStochasticEnd